Reference documentation for deal.II version 8.4.2
mg_transfer_matrix_free.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2016 - 2016 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/parallel_vector.h>
21 #include <deal.II/grid/tria.h>
22 #include <deal.II/grid/tria_iterator.h>
23 #include <deal.II/dofs/dof_tools.h>
24 #include <deal.II/fe/fe.h>
25 #include <deal.II/fe/fe_tools.h>
26 #include <deal.II/dofs/dof_accessor.h>
27 #include <deal.II/multigrid/mg_tools.h>
28 #include <deal.II/multigrid/mg_transfer_matrix_free.h>
29 
30 #include <deal.II/matrix_free/shape_info.h>
31 #include <deal.II/matrix_free/fe_evaluation.h>
32 
33 #include <algorithm>
34 
35 DEAL_II_NAMESPACE_OPEN
36 
37 
38 template<int dim, typename Number>
40  :
41  fe_degree(0),
42  element_is_continuous(false),
43  n_components(0),
44  n_child_cell_dofs(0)
45 {}
46 
47 
48 
49 template<int dim, typename Number>
51  :
52  fe_degree(0),
53  element_is_continuous(false),
54  n_components(0),
56 {
57  this->mg_constrained_dofs = &mg_c;
58 }
59 
60 
61 
62 template <int dim, typename Number>
64 {}
65 
66 
67 
68 template <int dim, typename Number>
70 (const MGConstrainedDoFs &mg_c)
71 {
72  this->mg_constrained_dofs = &mg_c;
73 }
74 
75 
76 
77 template <int dim, typename Number>
79 {
81  fe_degree = 0;
82  element_is_continuous = false;
83  n_components = 0;
85  level_dof_indices.clear();
86  parent_child_connect.clear();
87  n_owned_level_cells.clear();
90  weights_on_refined.clear();
91 }
92 
93 
94 
95 namespace
96 {
97  // given the collection of child cells in lexicographic ordering as seen
98  // from the parent, compute the first index of the given child
99  template <int dim>
100  unsigned int
101  compute_shift_within_children(const unsigned int child,
102  const unsigned int fe_shift_1d,
103  const unsigned int fe_degree)
104  {
105  // we put the degrees of freedom of all child cells in
106  // lexicographic ordering
107  unsigned int c_tensor_index[dim];
108  unsigned int tmp = child;
109  for (unsigned int d=0; d<dim; ++d)
110  {
111  c_tensor_index[d] = tmp % 2;
112  tmp /= 2;
113  }
114  const unsigned int n_child_dofs_1d = fe_degree + 1 + fe_shift_1d;
115  unsigned int factor = 1;
116  unsigned int shift = fe_shift_1d * c_tensor_index[0];
117  for (unsigned int d=1; d<dim; ++d)
118  {
119  factor *= n_child_dofs_1d;
120  shift = shift + factor * fe_shift_1d * c_tensor_index[d];
121  }
122  return shift;
123  }
124 
125 
126 
127  // puts the indices on the given child cell in lexicographic ordering with
128  // respect to the collection of all child cells as seen from the parent
129  template <int dim>
130  void add_child_indices(const unsigned int child,
131  const unsigned int fe_shift_1d,
132  const unsigned int fe_degree,
133  const std::vector<unsigned int> &lexicographic_numbering,
134  const std::vector<types::global_dof_index> &local_dof_indices,
135  types::global_dof_index *target_indices)
136  {
137  const unsigned int n_child_dofs_1d = fe_degree + 1 + fe_shift_1d;
138  const unsigned int shift =
139  compute_shift_within_children<dim>(child, fe_shift_1d, fe_degree);
140  const unsigned int n_components =
141  local_dof_indices.size()/Utilities::fixed_power<dim>(fe_degree+1);
142  types::global_dof_index *indices = target_indices + shift;
143  const unsigned int n_scalar_cell_dofs = Utilities::fixed_power<dim>(n_child_dofs_1d);
144  for (unsigned int c=0, m=0; c<n_components; ++c)
145  for (unsigned int k=0; k<(dim>2 ? (fe_degree+1) : 1); ++k)
146  for (unsigned int j=0; j<(dim>1 ? (fe_degree+1) : 1); ++j)
147  for (unsigned int i=0; i<(fe_degree+1); ++i, ++m)
148  {
149  const unsigned int index = c*n_scalar_cell_dofs+k*n_child_dofs_1d*
150  n_child_dofs_1d+j*n_child_dofs_1d+i;
151  Assert(indices[index] == numbers::invalid_dof_index ||
152  indices[index] == local_dof_indices[lexicographic_numbering[m]],
153  ExcInternalError());
154  indices[index] = local_dof_indices[lexicographic_numbering[m]];
155  }
156  }
157 
158 
159 
160  // initialize the vectors needed for the transfer (and merge with the
161  // content in copy_indices_global_mine)
162  template <typename Number>
163  void
164  reinit_ghosted_vector(const IndexSet &locally_owned,
165  std::vector<types::global_dof_index> &ghosted_level_dofs,
166  const MPI_Comm &communicator,
168  std::vector<std::pair<unsigned int,unsigned int> > &copy_indices_global_mine)
169  {
170  std::sort(ghosted_level_dofs.begin(), ghosted_level_dofs.end());
171  IndexSet ghosted_dofs(locally_owned.size());
172  ghosted_dofs.add_indices(ghosted_level_dofs.begin(),
173  std::unique(ghosted_level_dofs.begin(),
174  ghosted_level_dofs.end()));
175  ghosted_dofs.compress();
176 
177  // Add possible ghosts from the previous content in the vector
178  if (ghosted_level_vector.size() == locally_owned.size())
179  {
180  // shift the local number of the copy indices according to the new
181  // partitioner that we are going to use for the vector
182  const std_cxx11::shared_ptr<const Utilities::MPI::Partitioner> part
183  = ghosted_level_vector.get_partitioner();
184  ghosted_dofs.add_indices(part->ghost_indices());
185  for (unsigned int i=0; i<copy_indices_global_mine.size(); ++i)
186  copy_indices_global_mine[i].second =
187  locally_owned.n_elements() +
188  ghosted_dofs.index_within_set(part->local_to_global(copy_indices_global_mine[i].second));
189  }
190  ghosted_level_vector.reinit(locally_owned, ghosted_dofs, communicator);
191  }
192 
193  // Transform the ghost indices to local index space for the vector
194  void
195  copy_indices_to_mpi_local_numbers(const Utilities::MPI::Partitioner &part,
196  const std::vector<types::global_dof_index> &mine,
197  const std::vector<types::global_dof_index> &remote,
198  std::vector<unsigned int> &localized_indices)
199  {
200  localized_indices.resize(mine.size()+remote.size(),
202  for (unsigned int i=0; i<mine.size(); ++i)
203  if (mine[i] != numbers::invalid_dof_index)
204  localized_indices[i] = part.global_to_local(mine[i]);
205 
206  for (unsigned int i=0; i<remote.size(); ++i)
207  if (remote[i] != numbers::invalid_dof_index)
208  localized_indices[i+mine.size()] = part.global_to_local(remote[i]);
209  }
210 }
211 
212 
213 
214 template <int dim, typename Number>
216 (const DoFHandler<dim,dim> &mg_dof)
217 {
218  this->fill_and_communicate_copy_indices(mg_dof);
219 
220  // we collect all child DoFs of a mother cell together. For faster
221  // tensorized operations, we align the degrees of freedom
222  // lexicographically. We distinguish FE_Q elements and FE_DGQ elements
223 
224  const Triangulation<dim> &tria = mg_dof.get_triangulation();
225 
226  // ---------------------------- 1. Extract 1D info about the finite element
227  // step 1.1: create a 1D copy of the finite element from FETools where we
228  // substitute the template argument
229  AssertDimension(mg_dof.get_fe().n_base_elements(), 1);
230  std::string fe_name = mg_dof.get_fe().base_element(0).get_name();
231  {
232  const std::size_t template_starts = fe_name.find_first_of('<');
233  Assert (fe_name[template_starts+1] == (dim==1?'1':(dim==2?'2':'3')),
234  ExcInternalError());
235  fe_name[template_starts+1] = '1';
236  }
237  std_cxx11::shared_ptr<FiniteElement<1> > fe_1d
238  (FETools::get_fe_from_name<1>(fe_name));
239  const FiniteElement<1> &fe = *fe_1d;
240  unsigned int n_child_dofs_1d = numbers::invalid_unsigned_int;
241 
242  {
243  // currently, we have only FE_Q and FE_DGQ type elements implemented
244  n_components = mg_dof.get_fe().element_multiplicity(0);
245  AssertDimension(Utilities::fixed_power<dim>(fe.dofs_per_cell)*n_components,
246  mg_dof.get_fe().dofs_per_cell);
247  AssertDimension(fe.degree, mg_dof.get_fe().degree);
248  fe_degree = fe.degree;
250  Assert(fe.dofs_per_vertex < 2, ExcNotImplemented());
251 
252  // step 1.2: get renumbering of 1D basis functions to lexicographic
253  // numbers. The distinction according to fe.dofs_per_vertex is to support
254  // both continuous and discontinuous bases.
255  std::vector<unsigned int> renumbering(fe.dofs_per_cell);
256  {
258  renumbering[0] = 0;
259  for (unsigned int i=0; i<fe.dofs_per_line; ++i)
260  renumbering[i+fe.dofs_per_vertex] =
262  if (fe.dofs_per_vertex > 0)
263  renumbering[fe.dofs_per_cell-fe.dofs_per_vertex] = fe.dofs_per_vertex;
264  }
265 
266  // step 1.3: create a 1D quadrature formula from the finite element that
267  // collects the support points of the basis functions on the two children.
268  std::vector<Point<1> > basic_support_points = fe.get_unit_support_points();
269  Assert(fe.dofs_per_vertex == 0 || fe.dofs_per_vertex == 1,
270  ExcNotImplemented());
271  std::vector<Point<1> > points_refined(fe.dofs_per_vertex > 0 ?
272  (2 * fe.dofs_per_cell - 1) :
273  (2 * fe.dofs_per_cell));
274  const unsigned int shift = fe.dofs_per_cell - fe.dofs_per_vertex;
275  for (unsigned int c=0; c<GeometryInfo<1>::max_children_per_cell; ++c)
276  for (unsigned int j=0; j<basic_support_points.size(); ++j)
277  points_refined[shift*c+j][0] =
278  c*0.5 + 0.5 * basic_support_points[renumbering[j]][0];
279 
280  n_child_dofs_1d = points_refined.size();
281  n_child_cell_dofs = n_components*Utilities::fixed_power<dim>(n_child_dofs_1d);
282 
283  // step 1.4: evaluate the polynomials and store the data in ShapeInfo
284  const Quadrature<1> quadrature(points_refined);
285  shape_info.reinit(quadrature, mg_dof.get_fe(), 0);
286 
287  for (unsigned int c=0; c<GeometryInfo<1>::max_children_per_cell; ++c)
288  for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
289  for (unsigned int j=0; j<fe.dofs_per_cell; ++j)
290  Assert(std::abs(shape_info.shape_values[i*n_child_dofs_1d+j+c*shift][0] -
291  fe.get_prolongation_matrix(c)(renumbering[j],renumbering[i]))
292  < std::max(2.*(double)std::numeric_limits<Number>::epsilon(),1e-12),
293  ExcInternalError());
294  }
295 
296  // -------------- 2. Extract and match dof indices between child and parent
297  const unsigned int n_levels = tria.n_global_levels();
298  level_dof_indices.resize(n_levels);
299  parent_child_connect.resize(n_levels-1);
300  n_owned_level_cells.resize(n_levels-1);
301  std::vector<std::vector<unsigned int> > coarse_level_indices(n_levels-1);
302  for (unsigned int level=0; level<std::min(tria.n_levels(),n_levels-1); ++level)
303  coarse_level_indices[level].resize(tria.n_raw_cells(level),
305  std::vector<types::global_dof_index> local_dof_indices(mg_dof.get_fe().dofs_per_cell);
306  dirichlet_indices.resize(n_levels-1);
307 
308  // We use the vectors stored ghosted_level_vector in the base class for
309  // keeping ghosted transfer indices. To avoid keeping two very similar
310  // vectors, we merge them here.
311  if (this->ghosted_level_vector.max_level() != n_levels-1)
312  this->ghosted_level_vector.resize(0, n_levels-1);
313 
314  for (unsigned int level=n_levels-1; level > 0; --level)
315  {
316  unsigned int counter = 0;
317  std::vector<types::global_dof_index> global_level_dof_indices;
318  std::vector<types::global_dof_index> global_level_dof_indices_remote;
319  std::vector<types::global_dof_index> ghosted_level_dofs;
320  std::vector<types::global_dof_index> global_level_dof_indices_l0;
321  std::vector<types::global_dof_index> ghosted_level_dofs_l0;
322 
323  // step 2.1: loop over the cells on the coarse side
324  for (typename DoFHandler<dim>::cell_iterator cell = mg_dof.begin(level-1);
325  cell != mg_dof.end(level-1); ++cell)
326  {
327  // need to look into a cell if it has children and it is locally owned
328  if (!cell->has_children())
329  continue;
330 
331  bool consider_cell = false;
333  || cell->level_subdomain_id()==tria.locally_owned_subdomain()
334  )
335  consider_cell = true;
336 
337  // due to the particular way we store DoF indices (via children), we
338  // also need to add the DoF indices for coarse cells where we own at
339  // least one child
340  bool cell_is_remote = !consider_cell;
341  for (unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
342  if (cell->child(c)->level_subdomain_id()==tria.locally_owned_subdomain())
343  {
344  consider_cell = true;
345  break;
346  }
347 
348  if (!consider_cell)
349  continue;
350 
351  // step 2.2: loop through children and append the dof indices to the
352  // appropriate list. We need separate lists for the owned coarse
353  // cell case (which will be part of restriction/prolongation between
354  // level-1 and level) and the remote case (which needs to store DoF
355  // indices for the operations between level and level+1).
356  AssertDimension(cell->n_children(),
358  std::vector<types::global_dof_index> &next_indices =
359  cell_is_remote ? global_level_dof_indices_remote : global_level_dof_indices;
360  const std::size_t start_index = next_indices.size();
361  next_indices.resize(start_index + n_child_cell_dofs,
363  for (unsigned int c=0; c<GeometryInfo<dim>::max_children_per_cell; ++c)
364  {
365  if (cell_is_remote && cell->child(c)->level_subdomain_id() !=
367  continue;
368  cell->child(c)->get_mg_dof_indices(local_dof_indices);
369 
370  const IndexSet &owned_level_dofs = mg_dof.locally_owned_mg_dofs(level);
371  for (unsigned int i=0; i<local_dof_indices.size(); ++i)
372  if (!owned_level_dofs.is_element(local_dof_indices[i]))
373  ghosted_level_dofs.push_back(local_dof_indices[i]);
374 
375  add_child_indices<dim>(c, fe.dofs_per_cell - fe.dofs_per_vertex,
376  fe.degree, shape_info.lexicographic_numbering,
377  local_dof_indices,
378  &next_indices[start_index]);
379 
380  // step 2.3 store the connectivity to the parent
381  if (cell->child(c)->has_children() &&
383  || cell->child(c)->level_subdomain_id()==tria.locally_owned_subdomain()
384  ))
385  {
386  const unsigned int child_index = coarse_level_indices[level][cell->child(c)->index()];
387  AssertIndexRange(child_index, parent_child_connect[level].size());
388  unsigned int parent_index = counter;
389  // remote cells, i.e., cells where we work on a further
390  // level but are not treated on the current level, need to
391  // be placed at the end of the list; however, we do not yet
392  // know the exact position in the array, so shift their
393  // parent index by the number of cells so we can set the
394  // correct number after the end of this loop
395  if (cell_is_remote)
396  parent_index = start_index/n_child_cell_dofs + tria.n_cells(level);
397  parent_child_connect[level][child_index] =
398  std::make_pair(parent_index, c);
399  AssertIndexRange(mg_dof.get_fe().dofs_per_cell,
400  static_cast<unsigned short>(-1));
401 
402  // set Dirichlet boundary conditions (as a list of
403  // constrained DoFs) for the child
404  if (this->mg_constrained_dofs != 0)
405  for (unsigned int i=0; i<mg_dof.get_fe().dofs_per_cell; ++i)
406  if (this->mg_constrained_dofs->is_boundary_index(level, local_dof_indices[shape_info.lexicographic_numbering[i]]))
407  dirichlet_indices[level][child_index].push_back(i);
408  }
409  }
410  if (!cell_is_remote)
411  {
412  AssertIndexRange(static_cast<unsigned int>(cell->index()),
413  coarse_level_indices[level-1].size());
414  coarse_level_indices[level-1][cell->index()] = counter++;
415  }
416 
417  // step 2.4: include indices for the coarsest cells. we still insert
418  // the indices as if they were from a child in order to use the same
419  // code (the coarsest level does not matter much in terms of memory,
420  // so we gain in code simplicity)
421  if (level == 1 && !cell_is_remote)
422  {
423  cell->get_mg_dof_indices(local_dof_indices);
424 
425  const IndexSet &owned_level_dofs_l0 = mg_dof.locally_owned_mg_dofs(0);
426  for (unsigned int i=0; i<local_dof_indices.size(); ++i)
427  if (!owned_level_dofs_l0.is_element(local_dof_indices[i]))
428  ghosted_level_dofs_l0.push_back(local_dof_indices[i]);
429 
430  const std::size_t start_index = global_level_dof_indices_l0.size();
431  global_level_dof_indices_l0.resize(start_index+n_child_cell_dofs,
433  add_child_indices<dim>(0, fe.dofs_per_cell - fe.dofs_per_vertex,
434  fe.degree, shape_info.lexicographic_numbering,
435  local_dof_indices,
436  &global_level_dof_indices_l0[start_index]);
437 
438  dirichlet_indices[0].push_back(std::vector<unsigned short>());
439  if (this->mg_constrained_dofs != 0)
440  for (unsigned int i=0; i<mg_dof.get_fe().dofs_per_cell; ++i)
441  if (this->mg_constrained_dofs->is_boundary_index(0, local_dof_indices[shape_info.lexicographic_numbering[i]]))
442  dirichlet_indices[0].back().push_back(i);
443  }
444  }
445 
446  // step 2.5: store information about the current level and prepare the
447  // Dirichlet indices and parent-child relationship for the next coarser
448  // level
449  AssertDimension(counter*n_child_cell_dofs, global_level_dof_indices.size());
450  n_owned_level_cells[level-1] = counter;
451  dirichlet_indices[level-1].resize(counter);
452  parent_child_connect[level-1].
453  resize(counter, std::make_pair(numbers::invalid_unsigned_int,
455 
456  // step 2.6: put the cells with remotely owned parent to the end of the
457  // list (these are needed for the transfer from level to level+1 but not
458  // for the transfer from level-1 to level).
459  if (level < n_levels-1)
460  for (std::vector<std::pair<unsigned int,unsigned int> >::iterator
461  i=parent_child_connect[level].begin(); i!=parent_child_connect[level].end(); ++i)
462  if (i->first >= tria.n_cells(level))
463  {
464  i->first -= tria.n_cells(level);
465  i->first += counter;
466  }
467 
468  // step 2.7: Initialize the ghosted vector
469  const parallel::Triangulation<dim,dim> *ptria =
470  (dynamic_cast<const parallel::Triangulation<dim,dim>*> (&tria));
471  const MPI_Comm communicator =
472  ptria != 0 ? ptria->get_communicator() : MPI_COMM_SELF;
473 
474  reinit_ghosted_vector(mg_dof.locally_owned_mg_dofs(level),
475  ghosted_level_dofs, communicator,
476  this->ghosted_level_vector[level],
477  this->copy_indices_global_mine[level]);
478 
479  copy_indices_to_mpi_local_numbers(*this->ghosted_level_vector[level].get_partitioner(),
480  global_level_dof_indices,
481  global_level_dof_indices_remote,
482  level_dof_indices[level]);
483 
484  // step 2.8: Initialize the ghosted vector for level 0
485  if (level == 1)
486  {
487  for (unsigned int i = 0; i<parent_child_connect[0].size(); ++i)
488  parent_child_connect[0][i] = std::make_pair(i, 0U);
489 
490  reinit_ghosted_vector(mg_dof.locally_owned_mg_dofs(0),
491  ghosted_level_dofs_l0, communicator,
492  this->ghosted_level_vector[0],
493  this->copy_indices_global_mine[0]);
494 
495  copy_indices_to_mpi_local_numbers(*this->ghosted_level_vector[0].get_partitioner(),
496  global_level_dof_indices_l0,
497  std::vector<types::global_dof_index>(),
498  level_dof_indices[0]);
499  }
500  }
501 
502  // ------------------------ 3. compute weights to make restriction additive
503  //
504  // get the valence of the individual components and compute the weights as
505  // the inverse of the valence
506  weights_on_refined.resize(n_levels);
507  for (unsigned int level = 1; level<n_levels; ++level)
508  {
509  this->ghosted_level_vector[level] = 0;
510  for (unsigned int c=0; c<n_owned_level_cells[level-1]; ++c)
511  for (unsigned int j=0; j<n_child_cell_dofs; ++j)
512  this->ghosted_level_vector[level].local_element(level_dof_indices[level][n_child_cell_dofs*c+j]) += Number(1.);
513  this->ghosted_level_vector[level].compress(VectorOperation::add);
514  this->ghosted_level_vector[level].update_ghost_values();
515 
516  const unsigned int vec_size = VectorizedArray<Number>::n_array_elements;
517  std::vector<unsigned int> degree_to_3 (n_child_dofs_1d);
518  degree_to_3[0] = 0;
519  for (unsigned int i=1; i<n_child_dofs_1d-1; ++i)
520  degree_to_3[i] = 1;
521  degree_to_3.back() = 2;
522 
523  // we only store 3^dim weights because all dofs on a line have the same
524  // valence, and all dofs on a quad have the same valence.
525  weights_on_refined[level].resize(((n_owned_level_cells[level-1]+vec_size-1)/vec_size)*Utilities::fixed_power<dim>(3));
526  for (unsigned int c=0; c<n_owned_level_cells[level-1]; ++c)
527  {
528  const unsigned int comp = c/vec_size;
529  const unsigned int v = c%vec_size;
530 
531  for (unsigned int k=0, m=0; k<(dim>2 ? n_child_dofs_1d : 1); ++k)
532  for (unsigned int j=0; j<(dim>1 ? n_child_dofs_1d : 1); ++j)
533  {
534  unsigned int shift = 9*degree_to_3[k] + 3*degree_to_3[j];
535  for (unsigned int i=0; i<n_child_dofs_1d; ++i, ++m)
536  weights_on_refined[level][comp*Utilities::fixed_power<dim>(3)+shift+degree_to_3[i]][v] = Number(1.)/
537  this->ghosted_level_vector[level].local_element(level_dof_indices[level][n_child_cell_dofs*c+m]);
538  }
539  }
540  }
541 
543 }
544 
545 
546 
547 template <int dim, typename Number>
549 ::prolongate (const unsigned int to_level,
551  const parallel::distributed::Vector<Number> &src) const
552 {
553  Assert ((to_level >= 1) && (to_level<=level_dof_indices.size()),
554  ExcIndexRange (to_level, 1, level_dof_indices.size()+1));
555 
556  AssertDimension(this->ghosted_level_vector[to_level].local_size(),
557  dst.local_size());
558  AssertDimension(this->ghosted_level_vector[to_level-1].local_size(),
559  src.local_size());
560 
561  this->ghosted_level_vector[to_level-1] = src;
562  this->ghosted_level_vector[to_level-1].update_ghost_values();
563  this->ghosted_level_vector[to_level] = 0.;
564 
565  // the implementation in do_prolongate_add is templated in the degree of the
566  // element (for efficiency reasons), so we need to find the appropriate
567  // kernel here...
568  if (fe_degree == 0)
569  do_prolongate_add<0>(to_level, this->ghosted_level_vector[to_level],
570  this->ghosted_level_vector[to_level-1]);
571  else if (fe_degree == 1)
572  do_prolongate_add<1>(to_level, this->ghosted_level_vector[to_level],
573  this->ghosted_level_vector[to_level-1]);
574  else if (fe_degree == 2)
575  do_prolongate_add<2>(to_level, this->ghosted_level_vector[to_level],
576  this->ghosted_level_vector[to_level-1]);
577  else if (fe_degree == 3)
578  do_prolongate_add<3>(to_level, this->ghosted_level_vector[to_level],
579  this->ghosted_level_vector[to_level-1]);
580  else if (fe_degree == 4)
581  do_prolongate_add<4>(to_level, this->ghosted_level_vector[to_level],
582  this->ghosted_level_vector[to_level-1]);
583  else if (fe_degree == 5)
584  do_prolongate_add<5>(to_level, this->ghosted_level_vector[to_level],
585  this->ghosted_level_vector[to_level-1]);
586  else if (fe_degree == 6)
587  do_prolongate_add<6>(to_level, this->ghosted_level_vector[to_level],
588  this->ghosted_level_vector[to_level-1]);
589  else if (fe_degree == 7)
590  do_prolongate_add<7>(to_level, this->ghosted_level_vector[to_level],
591  this->ghosted_level_vector[to_level-1]);
592  else if (fe_degree == 8)
593  do_prolongate_add<8>(to_level, this->ghosted_level_vector[to_level],
594  this->ghosted_level_vector[to_level-1]);
595  else if (fe_degree == 9)
596  do_prolongate_add<9>(to_level, this->ghosted_level_vector[to_level],
597  this->ghosted_level_vector[to_level-1]);
598  else if (fe_degree == 10)
599  do_prolongate_add<10>(to_level, this->ghosted_level_vector[to_level],
600  this->ghosted_level_vector[to_level-1]);
601  else
602  AssertThrow(false, ExcNotImplemented("Only degrees 0 up to 10 implemented."));
603 
604  this->ghosted_level_vector[to_level].compress(VectorOperation::add);
605  dst = this->ghosted_level_vector[to_level];
606 }
607 
608 
609 
610 template <int dim, typename Number>
612 ::restrict_and_add (const unsigned int from_level,
614  const parallel::distributed::Vector<Number> &src) const
615 {
616  Assert ((from_level >= 1) && (from_level<=level_dof_indices.size()),
617  ExcIndexRange (from_level, 1, level_dof_indices.size()+1));
618 
619  AssertDimension(this->ghosted_level_vector[from_level].local_size(),
620  src.local_size());
621  AssertDimension(this->ghosted_level_vector[from_level-1].local_size(),
622  dst.local_size());
623 
624  this->ghosted_level_vector[from_level] = src;
625  this->ghosted_level_vector[from_level].update_ghost_values();
626  this->ghosted_level_vector[from_level-1] = 0.;
627 
628  if (fe_degree == 0)
629  do_restrict_add<0>(from_level, this->ghosted_level_vector[from_level-1],
630  this->ghosted_level_vector[from_level]);
631  else if (fe_degree == 1)
632  do_restrict_add<1>(from_level, this->ghosted_level_vector[from_level-1],
633  this->ghosted_level_vector[from_level]);
634  else if (fe_degree == 2)
635  do_restrict_add<2>(from_level, this->ghosted_level_vector[from_level-1],
636  this->ghosted_level_vector[from_level]);
637  else if (fe_degree == 3)
638  do_restrict_add<3>(from_level, this->ghosted_level_vector[from_level-1],
639  this->ghosted_level_vector[from_level]);
640  else if (fe_degree == 4)
641  do_restrict_add<4>(from_level, this->ghosted_level_vector[from_level-1],
642  this->ghosted_level_vector[from_level]);
643  else if (fe_degree == 5)
644  do_restrict_add<5>(from_level, this->ghosted_level_vector[from_level-1],
645  this->ghosted_level_vector[from_level]);
646  else if (fe_degree == 6)
647  do_restrict_add<6>(from_level, this->ghosted_level_vector[from_level-1],
648  this->ghosted_level_vector[from_level]);
649  else if (fe_degree == 7)
650  do_restrict_add<7>(from_level, this->ghosted_level_vector[from_level-1],
651  this->ghosted_level_vector[from_level]);
652  else if (fe_degree == 8)
653  do_restrict_add<8>(from_level, this->ghosted_level_vector[from_level-1],
654  this->ghosted_level_vector[from_level]);
655  else if (fe_degree == 9)
656  do_restrict_add<9>(from_level, this->ghosted_level_vector[from_level-1],
657  this->ghosted_level_vector[from_level]);
658  else if (fe_degree == 10)
659  do_restrict_add<10>(from_level, this->ghosted_level_vector[from_level-1],
660  this->ghosted_level_vector[from_level]);
661  else
662  AssertThrow(false, ExcNotImplemented("Only degrees 0 up to 10 implemented."));
663 
664  this->ghosted_level_vector[from_level-1].compress(VectorOperation::add);
665  dst += this->ghosted_level_vector[from_level-1];
666 }
667 
668 
669 
670 namespace
671 {
672  template <int dim, typename Eval, typename Number, bool prolongate>
673  void
674  perform_tensorized_op(const Eval &evaluator,
675  const unsigned int n_child_cell_dofs,
676  const unsigned int n_components,
678  {
679  AssertDimension(n_components * Eval::n_q_points, n_child_cell_dofs);
683 
684  for (unsigned int c=0; c<n_components; ++c)
685  {
686  // for the prolongate case, we go from dofs (living on the parent cell) to
687  // quads (living on all children) in the FEEvaluation terminology
688  if (dim == 1)
689  evaluator.template values<0,prolongate,false>(t0, t2);
690  else if (dim == 2)
691  {
692  evaluator.template values<0,prolongate,false>(t0, t1);
693  evaluator.template values<1,prolongate,false>(t1, t2);
694  }
695  else if (dim == 3)
696  {
697  evaluator.template values<0,prolongate,false>(t0, t2);
698  evaluator.template values<1,prolongate,false>(t2, t1);
699  evaluator.template values<2,prolongate,false>(t1, t2);
700  }
701  else
702  Assert(false, ExcNotImplemented());
703  if (prolongate)
704  {
705  t0 += Eval::dofs_per_cell;
706  t2 += Eval::n_q_points;
707  }
708  else
709  {
710  t0 += Eval::n_q_points;
711  t2 += Eval::dofs_per_cell;
712  }
713  }
714  }
715 
716  template <int dim, int degree, typename Number>
717  void weight_dofs_on_child (const VectorizedArray<Number> *weights,
718  const unsigned int n_components,
720  {
721  Assert(degree > 0, ExcNotImplemented());
722  const int loop_length = 2*degree+1;
723  unsigned int degree_to_3 [loop_length];
724  degree_to_3[0] = 0;
725  for (int i=1; i<loop_length-1; ++i)
726  degree_to_3[i] = 1;
727  degree_to_3[loop_length-1] = 2;
728  for (unsigned int c=0; c<n_components; ++c)
729  for (int k=0; k<(dim>2 ? loop_length : 1); ++k)
730  for (int j=0; j<(dim>1 ? loop_length : 1); ++j)
731  {
732  const unsigned int shift = 9*degree_to_3[k] + 3*degree_to_3[j];
733  data[0] *= weights[shift];
734  // loop bound as int avoids compiler warnings in case loop_length
735  // == 1 (polynomial degree 0)
736  for (int i=1; i<loop_length-1; ++i)
737  data[i] *= weights[shift+1];
738  data[loop_length-1] *= weights[shift+2];
739  data += loop_length;
740  }
741  }
742 }
743 
744 
745 
746 template <int dim, typename Number>
747 template <int degree>
749 ::do_prolongate_add (const unsigned int to_level,
751  const parallel::distributed::Vector<Number> &src) const
752 {
753  const unsigned int vec_size = VectorizedArray<Number>::n_array_elements;
754  const unsigned int n_child_dofs_1d = 2*(fe_degree+1) - element_is_continuous;
755  const unsigned int n_scalar_cell_dofs = Utilities::fixed_power<dim>(n_child_dofs_1d);
756  const unsigned int three_to_dim = Utilities::fixed_int_power<3,dim>::value;
757 
758  for (unsigned int cell=0; cell < n_owned_level_cells[to_level-1];
759  cell += vec_size)
760  {
761  const unsigned int n_chunks = cell+vec_size > n_owned_level_cells[to_level-1] ?
762  n_owned_level_cells[to_level-1] - cell : vec_size;
763 
764  // read from source vector
765  for (unsigned int v=0; v<n_chunks; ++v)
766  {
767  const unsigned int shift = compute_shift_within_children<dim>
768  (parent_child_connect[to_level-1][cell+v].second,
769  degree+1-element_is_continuous, degree);
770  const unsigned int *indices = &level_dof_indices[to_level-1][parent_child_connect[to_level-1][cell+v].first*n_child_cell_dofs+shift];
771  for (unsigned int c=0, m=0; c<n_components; ++c)
772  {
773  for (unsigned int k=0; k<(dim>2 ? (degree+1) : 1); ++k)
774  for (unsigned int j=0; j<(dim>1 ? (degree+1) : 1); ++j)
775  for (unsigned int i=0; i<(degree+1); ++i, ++m)
776  evaluation_data[m][v] =
777  src.local_element(indices[c*n_scalar_cell_dofs +
778  k*n_child_dofs_1d*n_child_dofs_1d+
779  j*n_child_dofs_1d+i]);
780 
781  // apply Dirichlet boundary conditions on parent cell
782  for (std::vector<unsigned short>::const_iterator i=dirichlet_indices[to_level-1][cell+v].begin(); i!=dirichlet_indices[to_level-1][cell+v].end(); ++i)
783  evaluation_data[*i][v] = 0.;
784  }
785  }
786 
787  // perform tensorized operation
788  Assert(shape_info.element_type ==
789  internal::MatrixFreeFunctions::tensor_symmetric, ExcNotImplemented());
791  {
792  AssertDimension(shape_info.shape_val_evenodd.size(),
793  (degree+1)*(degree+1));
794  typedef internal::EvaluatorTensorProduct<internal::evaluate_evenodd,dim,degree,2*degree+1,VectorizedArray<Number> > Evaluator;
795  Evaluator evaluator(shape_info.shape_val_evenodd,
796  shape_info.shape_val_evenodd,
797  shape_info.shape_val_evenodd);
798  perform_tensorized_op<dim,Evaluator,Number,true>(evaluator,
800  n_components,
802  weight_dofs_on_child<dim,degree,Number>(&weights_on_refined[to_level][(cell/vec_size)*three_to_dim],
803  n_components,
804  &evaluation_data[2*n_child_cell_dofs]);
805  }
806  else
807  {
808  AssertDimension(shape_info.shape_val_evenodd.size(),
809  (degree+1)*(degree+1));
810  typedef internal::EvaluatorTensorProduct<internal::evaluate_evenodd,dim,degree,2*degree+2,VectorizedArray<Number> > Evaluator;
811  Evaluator evaluator(shape_info.shape_val_evenodd,
812  shape_info.shape_val_evenodd,
813  shape_info.shape_val_evenodd);
814  perform_tensorized_op<dim,Evaluator,Number,true>(evaluator,
816  n_components,
818  }
819 
820  // write into dst vector
821  const unsigned int *indices = &level_dof_indices[to_level][cell*
823  for (unsigned int v=0; v<n_chunks; ++v)
824  {
825  for (unsigned int i=0; i<n_child_cell_dofs; ++i)
826  dst.local_element(indices[i]) += evaluation_data[2*n_child_cell_dofs+i][v];
827  indices += n_child_cell_dofs;
828  }
829  }
830 }
831 
832 
833 
834 template <int dim, typename Number>
835 template <int degree>
837 ::do_restrict_add (const unsigned int from_level,
839  const parallel::distributed::Vector<Number> &src) const
840 {
841  const unsigned int vec_size = VectorizedArray<Number>::n_array_elements;
842  const unsigned int n_child_dofs_1d = 2*(fe_degree+1) - element_is_continuous;
843  const unsigned int n_scalar_cell_dofs = Utilities::fixed_power<dim>(n_child_dofs_1d);
844  const unsigned int three_to_dim = Utilities::fixed_int_power<3,dim>::value;
845 
846  for (unsigned int cell=0; cell < n_owned_level_cells[from_level-1];
847  cell += vec_size)
848  {
849  const unsigned int n_chunks = cell+vec_size > n_owned_level_cells[from_level-1] ?
850  n_owned_level_cells[from_level-1] - cell : vec_size;
851 
852  // read from source vector
853  {
854  const unsigned int *indices = &level_dof_indices[from_level][cell*
856  for (unsigned int v=0; v<n_chunks; ++v)
857  {
858  for (unsigned int i=0; i<n_child_cell_dofs; ++i)
859  evaluation_data[i][v] = src.local_element(indices[i]);
860  indices += n_child_cell_dofs;
861  }
862  }
863 
864  // perform tensorized operation
865  Assert(shape_info.element_type ==
866  internal::MatrixFreeFunctions::tensor_symmetric, ExcNotImplemented());
868  {
869  AssertDimension(shape_info.shape_val_evenodd.size(),
870  (degree+1)*(degree+1));
871  typedef internal::EvaluatorTensorProduct<internal::evaluate_evenodd,dim,degree,2*degree+1,VectorizedArray<Number> > Evaluator;
872  Evaluator evaluator(shape_info.shape_val_evenodd,
873  shape_info.shape_val_evenodd,
874  shape_info.shape_val_evenodd);
875  weight_dofs_on_child<dim,degree,Number>(&weights_on_refined[from_level][(cell/vec_size)*three_to_dim],
876  n_components,
877  &evaluation_data[0]);
878  perform_tensorized_op<dim,Evaluator,Number,false>(evaluator,
880  n_components,
882  }
883  else
884  {
885  AssertDimension(shape_info.shape_val_evenodd.size(),
886  (degree+1)*(degree+1));
887  typedef internal::EvaluatorTensorProduct<internal::evaluate_evenodd,dim,degree,2*degree+2,VectorizedArray<Number> > Evaluator;
888  Evaluator evaluator(shape_info.shape_val_evenodd,
889  shape_info.shape_val_evenodd,
890  shape_info.shape_val_evenodd);
891  perform_tensorized_op<dim,Evaluator,Number,false>(evaluator,
893  n_components,
895  }
896 
897  // write into dst vector
898  for (unsigned int v=0; v<n_chunks; ++v)
899  {
900  const unsigned int shift = compute_shift_within_children<dim>
901  (parent_child_connect[from_level-1][cell+v].second,
902  degree+1-element_is_continuous, degree);
903  AssertIndexRange(parent_child_connect[from_level-1][cell+v].first*
905  level_dof_indices[from_level-1].size());
906  const unsigned int *indices = &level_dof_indices[from_level-1][parent_child_connect[from_level-1][cell+v].first*n_child_cell_dofs+shift];
907  for (unsigned int c=0, m=0; c<n_components; ++c)
908  {
909  // apply Dirichlet boundary conditions on parent cell
910  for (std::vector<unsigned short>::const_iterator i=dirichlet_indices[from_level-1][cell+v].begin(); i!=dirichlet_indices[from_level-1][cell+v].end(); ++i)
911  evaluation_data[2*n_child_cell_dofs+(*i)][v] = 0.;
912 
913  for (unsigned int k=0; k<(dim>2 ? (degree+1) : 1); ++k)
914  for (unsigned int j=0; j<(dim>1 ? (degree+1) : 1); ++j)
915  for (unsigned int i=0; i<(degree+1); ++i, ++m)
916  dst.local_element(indices[c*n_scalar_cell_dofs +
917  k*n_child_dofs_1d*n_child_dofs_1d+
918  j*n_child_dofs_1d+i])
920  }
921  }
922  }
923 }
924 
925 
926 
927 template <int dim, typename Number>
928 std::size_t
930 {
935  memory += shape_info.memory_consumption();
939  return memory;
940 }
941 
942 
943 // explicit instantiation
944 #include "mg_transfer_matrix_free.inst"
945 
946 
947 DEAL_II_NAMESPACE_CLOSE
unsigned int max_level() const
static const unsigned int invalid_unsigned_int
Definition: types.h:164
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1052
const types::subdomain_id invalid_subdomain_id
Definition: types.h:239
std::vector< unsigned int > n_owned_level_cells
cell_iterator begin(const unsigned int level=0) const
Definition: dof_handler.cc:833
AlignedVector< VectorizedArray< Number > > evaluation_data
unsigned int n_raw_cells(const unsigned int level) const
Definition: tria.cc:11039
std::vector< std::vector< std::pair< unsigned int, unsigned int > > > copy_indices_global_mine
Definition: mg_transfer.h:372
std_cxx11::shared_ptr< const Utilities::MPI::Partitioner > get_partitioner() const
unsigned int n_cells() const
Definition: tria.cc:10966
std::size_t memory_consumption() const
bool is_boundary_index(const unsigned int level, const types::global_dof_index index) const
DEAL_VOLATILE unsigned int counter
Definition: subscriptor.h:184
std::vector< std::vector< std::vector< unsigned short > > > dirichlet_indices
cell_iterator end() const
Definition: dof_handler.cc:861
#define AssertIndexRange(index, range)
Definition: exceptions.h:1081
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
Definition: index_set.h:1291
const unsigned int degree
Definition: fe_base.h:299
#define AssertThrow(cond, exc)
Definition: exceptions.h:358
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe.cc:326
unsigned int global_to_local(const types::global_dof_index global_index) const
const FiniteElement< dim, spacedim > & get_fe() const
const unsigned int dofs_per_line
Definition: fe_base.h:232
SmartPointer< const MGConstrainedDoFs, MGLevelGlobalTransfer< parallel::distributed::Vector< Number > > > mg_constrained_dofs
Definition: mg_transfer.h:401
unsigned int n_levels() const
virtual void restrict_and_add(const unsigned int from_level, parallel::distributed::Vector< Number > &dst, const parallel::distributed::Vector< Number > &src) const
size_type size() const
Definition: index_set.h:1247
Number local_element(const size_type local_index) const
void fill_and_communicate_copy_indices(const DoFHandler< dim, spacedim > &mg_dof)
std::vector< std::vector< std::pair< unsigned int, unsigned int > > > parent_child_connect
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:294
virtual MPI_Comm get_communicator() const
Definition: tria_base.cc:135
std::vector< AlignedVector< VectorizedArray< Number > > > weights_on_refined
void resize(const size_type size_in, const T &init=T())
MGLevelObject< parallel::distributed::Vector< Number > > ghosted_level_vector
Definition: mg_transfer.h:414
const IndexSet & locally_owned_mg_dofs(const unsigned int level) const
std::vector< std::vector< unsigned int > > level_dof_indices
const unsigned int dofs_per_cell
Definition: fe_base.h:283
void reinit(const size_type size, const bool omit_zeroing_entries=false)
const std::vector< Point< dim > > & get_unit_support_points() const
Definition: fe.cc:956
std_cxx11::enable_if< std_cxx11::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
virtual void prolongate(const unsigned int to_level, parallel::distributed::Vector< Number > &dst, const parallel::distributed::Vector< Number > &src) const
virtual unsigned int n_global_levels() const
void do_prolongate_add(const unsigned int to_level, parallel::distributed::Vector< Number > &dst, const parallel::distributed::Vector< Number > &src) const
const Triangulation< dim, spacedim > & get_triangulation() const
internal::MatrixFreeFunctions::ShapeInfo< Number > shape_info
const unsigned int dofs_per_vertex
Definition: fe_base.h:226
bool is_element(const size_type index) const
Definition: index_set.h:1317
void resize(const unsigned int new_minlevel, const unsigned int new_maxlevel)
const types::global_dof_index invalid_dof_index
Definition: types.h:178
void do_restrict_add(const unsigned int from_level, parallel::distributed::Vector< Number > &dst, const parallel::distributed::Vector< Number > &src) const
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
size_type local_size() const
size_type n_elements() const
Definition: index_set.h:1380
void build(const DoFHandler< dim, dim > &mg_dof)
virtual types::subdomain_id locally_owned_subdomain() const
Definition: tria.cc:11575