Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
mg_transfer_matrix_free.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2016 - 2020 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #include <deal.II/base/function.h>
18 #include <deal.II/base/logstream.h>
20 
22 #include <deal.II/dofs/dof_tools.h>
23 
24 #include <deal.II/fe/fe.h>
25 
26 #include <deal.II/grid/tria.h>
28 
30 
32 
36 
37 #include <algorithm>
38 
40 
41 
42 template <int dim, typename Number>
44  : fe_degree(0)
45  , element_is_continuous(false)
46  , n_components(0)
47  , n_child_cell_dofs(0)
48 {}
49 
50 
51 
52 template <int dim, typename Number>
54  const MGConstrainedDoFs &mg_c)
55  : fe_degree(0)
56  , element_is_continuous(false)
57  , n_components(0)
58  , n_child_cell_dofs(0)
59 {
60  this->mg_constrained_dofs = &mg_c;
61 }
62 
63 
64 
65 template <int dim, typename Number>
66 void
68  const MGConstrainedDoFs &mg_c)
69 {
70  this->mg_constrained_dofs = &mg_c;
71 }
72 
73 
74 
75 template <int dim, typename Number>
76 void
78 {
81  fe_degree = 0;
82  element_is_continuous = false;
83  n_components = 0;
84  n_child_cell_dofs = 0;
85  level_dof_indices.clear();
86  parent_child_connect.clear();
87  dirichlet_indices.clear();
88  n_owned_level_cells.clear();
89  prolongation_matrix_1d.clear();
90  evaluation_data.clear();
91  weights_on_refined.clear();
92 }
93 
94 
95 template <int dim, typename Number>
96 void
98  const DoFHandler<dim, dim> &dof_handler,
99  const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
100  &external_partitioners)
101 {
102  this->fill_and_communicate_copy_indices(dof_handler);
103 
104  vector_partitioners.resize(0,
105  dof_handler.get_triangulation().n_global_levels() -
106  1);
107  for (unsigned int level = 0; level <= this->ghosted_level_vector.max_level();
108  ++level)
109  vector_partitioners[level] =
110  this->ghosted_level_vector[level].get_partitioner();
111 
112  std::vector<std::vector<Number>> weights_unvectorized;
113 
115 
116  internal::MGTransfer::setup_transfer<dim, Number>(
117  dof_handler,
118  this->mg_constrained_dofs,
119  external_partitioners,
120  elem_info,
121  level_dof_indices,
122  parent_child_connect,
123  n_owned_level_cells,
124  dirichlet_indices,
125  weights_unvectorized,
126  this->copy_indices_global_mine,
127  vector_partitioners);
128 
129  // reinit the ghosted level vector in case its partitioner differs from the
130  // one the user intends to use externally
131  this->ghosted_level_vector.resize(0, vector_partitioners.max_level());
132  for (unsigned int level = 0; level <= vector_partitioners.max_level();
133  ++level)
134  if (external_partitioners.size() == vector_partitioners.max_level() + 1 &&
135  external_partitioners[level].get() == vector_partitioners[level].get())
136  this->ghosted_level_vector[level].reinit(0);
137  else
138  this->ghosted_level_vector[level].reinit(vector_partitioners[level]);
139 
140  // unpack element info data
141  fe_degree = elem_info.fe_degree;
142  element_is_continuous = elem_info.element_is_continuous;
143  n_components = elem_info.n_components;
144  n_child_cell_dofs = elem_info.n_child_cell_dofs;
145 
146  // duplicate and put into vectorized array
147  prolongation_matrix_1d.resize(elem_info.prolongation_matrix_1d.size());
148  for (unsigned int i = 0; i < elem_info.prolongation_matrix_1d.size(); i++)
149  prolongation_matrix_1d[i] = elem_info.prolongation_matrix_1d[i];
150 
151  // reshuffle into aligned vector of vectorized arrays
152  const unsigned int vec_size = VectorizedArray<Number>::size();
153  const unsigned int n_levels =
154  dof_handler.get_triangulation().n_global_levels();
155 
156  const unsigned int n_weights_per_cell = Utilities::fixed_power<dim>(3);
157  weights_on_refined.resize(n_levels - 1);
158  for (unsigned int level = 1; level < n_levels; ++level)
159  {
160  weights_on_refined[level - 1].resize(
161  ((n_owned_level_cells[level - 1] + vec_size - 1) / vec_size) *
162  n_weights_per_cell);
163 
164  for (unsigned int c = 0; c < n_owned_level_cells[level - 1]; ++c)
165  {
166  const unsigned int comp = c / vec_size;
167  const unsigned int v = c % vec_size;
168  for (unsigned int i = 0; i < n_weights_per_cell; ++i)
169  {
170  weights_on_refined[level - 1][comp * n_weights_per_cell + i][v] =
171  weights_unvectorized[level - 1][c * n_weights_per_cell + i];
172  }
173  }
174  }
175 
176  evaluation_data.resize(n_child_cell_dofs);
177 }
178 
179 
180 
181 template <int dim, typename Number>
182 void
184  const unsigned int to_level,
187 {
188  Assert((to_level >= 1) && (to_level <= level_dof_indices.size()),
189  ExcIndexRange(to_level, 1, level_dof_indices.size() + 1));
190 
191  const bool src_inplace = src.get_partitioner().get() ==
192  this->vector_partitioners[to_level - 1].get();
193  if (src_inplace == false)
194  {
195  if (this->ghosted_level_vector[to_level - 1].get_partitioner().get() !=
196  this->vector_partitioners[to_level - 1].get())
197  this->ghosted_level_vector[to_level - 1].reinit(
198  this->vector_partitioners[to_level - 1]);
199  this->ghosted_level_vector[to_level - 1].copy_locally_owned_data_from(
200  src);
201  }
202 
203  const bool dst_inplace =
204  dst.get_partitioner().get() == this->vector_partitioners[to_level].get();
205  if (dst_inplace == false)
206  {
207  if (this->ghosted_level_vector[to_level].get_partitioner().get() !=
208  this->vector_partitioners[to_level].get())
209  this->ghosted_level_vector[to_level].reinit(
210  this->vector_partitioners[to_level]);
211  AssertDimension(this->ghosted_level_vector[to_level].local_size(),
212  dst.local_size());
213  this->ghosted_level_vector[to_level] = 0.;
214  }
215  else
216  dst = 0;
217 
219  src_inplace ? src : this->ghosted_level_vector[to_level - 1];
221  dst_inplace ? dst : this->ghosted_level_vector[to_level];
222 
223  src_vec.update_ghost_values();
224  // the implementation in do_prolongate_add is templated in the degree of the
225  // element (for efficiency reasons), so we need to find the appropriate
226  // kernel here...
227  if (fe_degree == 0)
228  do_prolongate_add<0>(to_level, dst_vec, src_vec);
229  else if (fe_degree == 1)
230  do_prolongate_add<1>(to_level, dst_vec, src_vec);
231  else if (fe_degree == 2)
232  do_prolongate_add<2>(to_level, dst_vec, src_vec);
233  else if (fe_degree == 3)
234  do_prolongate_add<3>(to_level, dst_vec, src_vec);
235  else if (fe_degree == 4)
236  do_prolongate_add<4>(to_level, dst_vec, src_vec);
237  else if (fe_degree == 5)
238  do_prolongate_add<5>(to_level, dst_vec, src_vec);
239  else if (fe_degree == 6)
240  do_prolongate_add<6>(to_level, dst_vec, src_vec);
241  else if (fe_degree == 7)
242  do_prolongate_add<7>(to_level, dst_vec, src_vec);
243  else if (fe_degree == 8)
244  do_prolongate_add<8>(to_level, dst_vec, src_vec);
245  else if (fe_degree == 9)
246  do_prolongate_add<9>(to_level, dst_vec, src_vec);
247  else if (fe_degree == 10)
248  do_prolongate_add<10>(to_level, dst_vec, src_vec);
249  else
250  do_prolongate_add<-1>(to_level, dst_vec, src_vec);
251 
252  dst_vec.compress(VectorOperation::add);
253  if (dst_inplace == false)
254  dst = dst_vec;
255 
256  if (src_inplace == true)
257  src.zero_out_ghosts();
258 }
259 
260 
261 
262 template <int dim, typename Number>
263 void
265  const unsigned int from_level,
268 {
269  Assert((from_level >= 1) && (from_level <= level_dof_indices.size()),
270  ExcIndexRange(from_level, 1, level_dof_indices.size() + 1));
271 
272  const bool src_inplace =
273  src.get_partitioner().get() == this->vector_partitioners[from_level].get();
274  if (src_inplace == false)
275  {
276  if (this->ghosted_level_vector[from_level].get_partitioner().get() !=
277  this->vector_partitioners[from_level].get())
278  this->ghosted_level_vector[from_level].reinit(
279  this->vector_partitioners[from_level]);
280  this->ghosted_level_vector[from_level].copy_locally_owned_data_from(src);
281  }
282 
283  const bool dst_inplace = dst.get_partitioner().get() ==
284  this->vector_partitioners[from_level - 1].get();
285  if (dst_inplace == false)
286  {
287  if (this->ghosted_level_vector[from_level - 1].get_partitioner().get() !=
288  this->vector_partitioners[from_level - 1].get())
289  this->ghosted_level_vector[from_level - 1].reinit(
290  this->vector_partitioners[from_level - 1]);
291  AssertDimension(this->ghosted_level_vector[from_level - 1].local_size(),
292  dst.local_size());
293  this->ghosted_level_vector[from_level - 1] = 0.;
294  }
295 
297  src_inplace ? src : this->ghosted_level_vector[from_level];
299  dst_inplace ? dst : this->ghosted_level_vector[from_level - 1];
300 
301  src_vec.update_ghost_values();
302 
303  if (fe_degree == 0)
304  do_restrict_add<0>(from_level, dst_vec, src_vec);
305  else if (fe_degree == 1)
306  do_restrict_add<1>(from_level, dst_vec, src_vec);
307  else if (fe_degree == 2)
308  do_restrict_add<2>(from_level, dst_vec, src_vec);
309  else if (fe_degree == 3)
310  do_restrict_add<3>(from_level, dst_vec, src_vec);
311  else if (fe_degree == 4)
312  do_restrict_add<4>(from_level, dst_vec, src_vec);
313  else if (fe_degree == 5)
314  do_restrict_add<5>(from_level, dst_vec, src_vec);
315  else if (fe_degree == 6)
316  do_restrict_add<6>(from_level, dst_vec, src_vec);
317  else if (fe_degree == 7)
318  do_restrict_add<7>(from_level, dst_vec, src_vec);
319  else if (fe_degree == 8)
320  do_restrict_add<8>(from_level, dst_vec, src_vec);
321  else if (fe_degree == 9)
322  do_restrict_add<9>(from_level, dst_vec, src_vec);
323  else if (fe_degree == 10)
324  do_restrict_add<10>(from_level, dst_vec, src_vec);
325  else
326  // go to the non-templated version of the evaluator
327  do_restrict_add<-1>(from_level, dst_vec, src_vec);
328 
329  dst_vec.compress(VectorOperation::add);
330  if (dst_inplace == false)
331  dst += dst_vec;
332 
333  if (src_inplace == true)
334  src.zero_out_ghosts();
335 }
336 
337 
338 
339 namespace
340 {
341  template <int dim, int degree, typename Number>
342  void
343  weight_dofs_on_child(const VectorizedArray<Number> *weights,
344  const unsigned int n_components,
345  const unsigned int fe_degree,
347  {
348  Assert(fe_degree > 0, ExcNotImplemented());
349  Assert(fe_degree < 100, ExcNotImplemented());
350  const int loop_length = degree != -1 ? 2 * degree + 1 : 2 * fe_degree + 1;
351  unsigned int degree_to_3[100];
352  degree_to_3[0] = 0;
353  for (int i = 1; i < loop_length - 1; ++i)
354  degree_to_3[i] = 1;
355  degree_to_3[loop_length - 1] = 2;
356  for (unsigned int c = 0; c < n_components; ++c)
357  for (int k = 0; k < (dim > 2 ? loop_length : 1); ++k)
358  for (int j = 0; j < (dim > 1 ? loop_length : 1); ++j)
359  {
360  const unsigned int shift = 9 * degree_to_3[k] + 3 * degree_to_3[j];
361  data[0] *= weights[shift];
362  // loop bound as int avoids compiler warnings in case loop_length
363  // == 1 (polynomial degree 0)
364  for (int i = 1; i < loop_length - 1; ++i)
365  data[i] *= weights[shift + 1];
366  data[loop_length - 1] *= weights[shift + 2];
367  data += loop_length;
368  }
369  }
370 } // namespace
371 
372 
373 
374 template <int dim, typename Number>
375 template <int degree>
376 void
378  const unsigned int to_level,
381 {
382  const unsigned int vec_size = VectorizedArray<Number>::size();
383  const unsigned int degree_size = (degree > -1 ? degree : fe_degree) + 1;
384  const unsigned int n_child_dofs_1d = 2 * degree_size - element_is_continuous;
385  const unsigned int n_scalar_cell_dofs =
386  Utilities::fixed_power<dim>(n_child_dofs_1d);
387  constexpr unsigned int three_to_dim = Utilities::pow(3, dim);
388 
389  // If we have user defined MG constraints, we must create
390  // a non-const, ghosted version of the source vector to distribute
391  // constraints.
392  const LinearAlgebra::distributed::Vector<Number> *to_use = &src;
394  if (this->mg_constrained_dofs != nullptr &&
395  this->mg_constrained_dofs->get_user_constraint_matrix(to_level - 1)
396  .get_local_lines()
397  .size() > 0)
398  {
400 
401  // Distribute any user defined constraints
402  this->mg_constrained_dofs->get_user_constraint_matrix(to_level - 1)
403  .distribute(copy_src);
404 
405  // Re-initialize new ghosted vector with correct constraints
406  new_src.reinit(copy_src);
407  new_src = copy_src;
408  new_src.update_ghost_values();
409 
410  to_use = &new_src;
411  }
412 
413  for (unsigned int cell = 0; cell < n_owned_level_cells[to_level - 1];
414  cell += vec_size)
415  {
416  const unsigned int n_lanes =
417  cell + vec_size > n_owned_level_cells[to_level - 1] ?
418  n_owned_level_cells[to_level - 1] - cell :
419  vec_size;
420 
421  // read from source vector
422  for (unsigned int v = 0; v < n_lanes; ++v)
423  {
424  const unsigned int shift =
425  internal::MGTransfer::compute_shift_within_children<dim>(
426  parent_child_connect[to_level - 1][cell + v].second,
427  fe_degree + 1 - element_is_continuous,
428  fe_degree);
429  const unsigned int *indices =
430  &level_dof_indices[to_level - 1]
431  [parent_child_connect[to_level - 1][cell + v]
432  .first *
433  n_child_cell_dofs +
434  shift];
435  for (unsigned int c = 0, m = 0; c < n_components; ++c)
436  {
437  for (unsigned int k = 0; k < (dim > 2 ? degree_size : 1); ++k)
438  for (unsigned int j = 0; j < (dim > 1 ? degree_size : 1); ++j)
439  for (unsigned int i = 0; i < degree_size; ++i, ++m)
440  evaluation_data[m][v] = to_use->local_element(
441  indices[c * n_scalar_cell_dofs +
442  k * n_child_dofs_1d * n_child_dofs_1d +
443  j * n_child_dofs_1d + i]);
444 
445  // apply Dirichlet boundary conditions on parent cell
446  for (std::vector<unsigned short>::const_iterator i =
447  dirichlet_indices[to_level - 1][cell + v].begin();
448  i != dirichlet_indices[to_level - 1][cell + v].end();
449  ++i)
450  evaluation_data[*i][v] = 0.;
451  }
452  }
453 
454  AssertDimension(prolongation_matrix_1d.size(),
455  degree_size * n_child_dofs_1d);
456  // perform tensorized operation
457  if (element_is_continuous)
458  {
459  // must go through the components backwards because we want to write
460  // the output to the same array as the input
461  for (int c = n_components - 1; c >= 0; --c)
463  dim,
464  degree + 1,
465  2 * degree + 1,
466  1,
469  do_forward(prolongation_matrix_1d,
470  evaluation_data.begin() +
471  c * Utilities::fixed_power<dim>(degree_size),
472  evaluation_data.begin() + c * n_scalar_cell_dofs,
473  fe_degree + 1,
474  2 * fe_degree + 1);
475  weight_dofs_on_child<dim, degree, Number>(
476  &weights_on_refined[to_level - 1][(cell / vec_size) * three_to_dim],
477  n_components,
478  fe_degree,
479  evaluation_data.begin());
480  }
481  else
482  {
483  for (int c = n_components - 1; c >= 0; --c)
485  dim,
486  degree + 1,
487  2 * degree + 2,
488  1,
491  do_forward(prolongation_matrix_1d,
492  evaluation_data.begin() +
493  c * Utilities::fixed_power<dim>(degree_size),
494  evaluation_data.begin() + c * n_scalar_cell_dofs,
495  fe_degree + 1,
496  2 * fe_degree + 2);
497  }
498 
499  // write into dst vector
500  const unsigned int *indices =
501  &level_dof_indices[to_level][cell * n_child_cell_dofs];
502  for (unsigned int v = 0; v < n_lanes; ++v)
503  {
504  for (unsigned int i = 0; i < n_child_cell_dofs; ++i)
505  dst.local_element(indices[i]) += evaluation_data[i][v];
506  indices += n_child_cell_dofs;
507  }
508  }
509 }
510 
511 
512 
513 template <int dim, typename Number>
514 template <int degree>
515 void
517  const unsigned int from_level,
520 {
521  const unsigned int vec_size = VectorizedArray<Number>::size();
522  const unsigned int degree_size = (degree > -1 ? degree : fe_degree) + 1;
523  const unsigned int n_child_dofs_1d = 2 * degree_size - element_is_continuous;
524  const unsigned int n_scalar_cell_dofs =
525  Utilities::fixed_power<dim>(n_child_dofs_1d);
526  constexpr unsigned int three_to_dim = Utilities::pow(3, dim);
527 
528  for (unsigned int cell = 0; cell < n_owned_level_cells[from_level - 1];
529  cell += vec_size)
530  {
531  const unsigned int n_lanes =
532  cell + vec_size > n_owned_level_cells[from_level - 1] ?
533  n_owned_level_cells[from_level - 1] - cell :
534  vec_size;
535 
536  // read from source vector
537  {
538  const unsigned int *indices =
539  &level_dof_indices[from_level][cell * n_child_cell_dofs];
540  for (unsigned int v = 0; v < n_lanes; ++v)
541  {
542  for (unsigned int i = 0; i < n_child_cell_dofs; ++i)
543  evaluation_data[i][v] = src.local_element(indices[i]);
544  indices += n_child_cell_dofs;
545  }
546  }
547 
548  AssertDimension(prolongation_matrix_1d.size(),
549  degree_size * n_child_dofs_1d);
550  // perform tensorized operation
551  if (element_is_continuous)
552  {
553  weight_dofs_on_child<dim, degree, Number>(
554  &weights_on_refined[from_level - 1]
555  [(cell / vec_size) * three_to_dim],
556  n_components,
557  fe_degree,
558  evaluation_data.data());
559  for (unsigned int c = 0; c < n_components; ++c)
561  dim,
562  degree + 1,
563  2 * degree + 1,
564  1,
567  do_backward(prolongation_matrix_1d,
568  false,
569  evaluation_data.begin() + c * n_scalar_cell_dofs,
570  evaluation_data.begin() +
571  c * Utilities::fixed_power<dim>(degree_size),
572  fe_degree + 1,
573  2 * fe_degree + 1);
574  }
575  else
576  {
577  for (unsigned int c = 0; c < n_components; ++c)
579  dim,
580  degree + 1,
581  2 * degree + 2,
582  1,
585  do_backward(prolongation_matrix_1d,
586  false,
587  evaluation_data.begin() + c * n_scalar_cell_dofs,
588  evaluation_data.begin() +
589  c * Utilities::fixed_power<dim>(degree_size),
590  fe_degree + 1,
591  2 * fe_degree + 2);
592  }
593 
594  // write into dst vector
595  for (unsigned int v = 0; v < n_lanes; ++v)
596  {
597  const unsigned int shift =
598  internal::MGTransfer::compute_shift_within_children<dim>(
599  parent_child_connect[from_level - 1][cell + v].second,
600  fe_degree + 1 - element_is_continuous,
601  fe_degree);
603  parent_child_connect[from_level - 1][cell + v].first *
604  n_child_cell_dofs +
605  n_child_cell_dofs - 1,
606  level_dof_indices[from_level - 1].size());
607  const unsigned int *indices =
608  &level_dof_indices[from_level - 1]
609  [parent_child_connect[from_level - 1][cell + v]
610  .first *
611  n_child_cell_dofs +
612  shift];
613  for (unsigned int c = 0, m = 0; c < n_components; ++c)
614  {
615  // apply Dirichlet boundary conditions on parent cell
616  for (std::vector<unsigned short>::const_iterator i =
617  dirichlet_indices[from_level - 1][cell + v].begin();
618  i != dirichlet_indices[from_level - 1][cell + v].end();
619  ++i)
620  evaluation_data[*i][v] = 0.;
621 
622  for (unsigned int k = 0; k < (dim > 2 ? degree_size : 1); ++k)
623  for (unsigned int j = 0; j < (dim > 1 ? degree_size : 1); ++j)
624  for (unsigned int i = 0; i < degree_size; ++i, ++m)
625  dst.local_element(
626  indices[c * n_scalar_cell_dofs +
627  k * n_child_dofs_1d * n_child_dofs_1d +
628  j * n_child_dofs_1d + i]) +=
629  evaluation_data[m][v];
630  }
631  }
632  }
633 }
634 
635 
636 
637 template <int dim, typename Number>
638 std::size_t
640 {
641  std::size_t memory = MGLevelGlobalTransfer<
643  memory += MemoryConsumption::memory_consumption(level_dof_indices);
644  memory += MemoryConsumption::memory_consumption(parent_child_connect);
645  memory += MemoryConsumption::memory_consumption(n_owned_level_cells);
646  memory += MemoryConsumption::memory_consumption(prolongation_matrix_1d);
647  memory += MemoryConsumption::memory_consumption(evaluation_data);
648  memory += MemoryConsumption::memory_consumption(weights_on_refined);
649  memory += MemoryConsumption::memory_consumption(dirichlet_indices);
650  return memory;
651 }
652 
653 
654 
655 template <int dim, typename Number>
657  const MGConstrainedDoFs &mg_c)
658  : same_for_all(true)
659 {
660  matrix_free_transfer_vector.emplace_back(mg_c);
661 }
662 
663 
664 
665 template <int dim, typename Number>
667  const std::vector<MGConstrainedDoFs> &mg_c)
668  : same_for_all(false)
669 {
670  for (const auto &constrained_block_dofs : mg_c)
671  matrix_free_transfer_vector.emplace_back(constrained_block_dofs);
672 }
673 
674 
675 
676 template <int dim, typename Number>
677 void
679  const MGConstrainedDoFs &mg_c)
680 {
681  Assert(same_for_all,
682  ExcMessage("This object was initialized with support for usage with "
683  "one DoFHandler for each block, but this method assumes "
684  "that the same DoFHandler is used for all the blocks!"));
685  AssertDimension(matrix_free_transfer_vector.size(), 1);
686 
687  matrix_free_transfer_vector[0].initialize_constraints(mg_c);
688 }
689 
690 
691 
692 template <int dim, typename Number>
693 void
695  const std::vector<MGConstrainedDoFs> &mg_c)
696 {
697  Assert(!same_for_all,
698  ExcMessage("This object was initialized with support for using "
699  "the same DoFHandler for all the blocks, but this "
700  "method assumes that there is a separate DoFHandler "
701  "for each block!"));
702  AssertDimension(matrix_free_transfer_vector.size(), mg_c.size());
703 
704  for (unsigned int i = 0; i < mg_c.size(); ++i)
705  matrix_free_transfer_vector[i].initialize_constraints(mg_c[i]);
706 }
707 
708 
709 
710 template <int dim, typename Number>
711 void
713 {
714  matrix_free_transfer_vector.clear();
715 }
716 
717 
718 
719 template <int dim, typename Number>
720 void
722  const DoFHandler<dim, dim> &dof_handler)
723 {
724  AssertDimension(matrix_free_transfer_vector.size(), 1);
725  matrix_free_transfer_vector[0].build(dof_handler);
726 }
727 
728 
729 
730 template <int dim, typename Number>
731 void
733  const std::vector<const DoFHandler<dim, dim> *> &dof_handler)
734 {
735  AssertDimension(matrix_free_transfer_vector.size(), dof_handler.size());
736  for (unsigned int i = 0; i < dof_handler.size(); ++i)
737  matrix_free_transfer_vector[i].build(*dof_handler[i]);
738 }
739 
740 
741 
742 template <int dim, typename Number>
743 void
745  const unsigned int to_level,
748 {
749  const unsigned int n_blocks = src.n_blocks();
751 
752  if (!same_for_all)
753  AssertDimension(matrix_free_transfer_vector.size(), n_blocks);
754 
755  for (unsigned int b = 0; b < n_blocks; ++b)
756  {
757  const unsigned int data_block = same_for_all ? 0 : b;
758  matrix_free_transfer_vector[data_block].prolongate(to_level,
759  dst.block(b),
760  src.block(b));
761  }
762 }
763 
764 
765 
766 template <int dim, typename Number>
767 void
769  const unsigned int from_level,
772 {
773  const unsigned int n_blocks = src.n_blocks();
775 
776  if (!same_for_all)
777  AssertDimension(matrix_free_transfer_vector.size(), n_blocks);
778 
779  for (unsigned int b = 0; b < n_blocks; ++b)
780  {
781  const unsigned int data_block = same_for_all ? 0 : b;
782  matrix_free_transfer_vector[data_block].restrict_and_add(from_level,
783  dst.block(b),
784  src.block(b));
785  }
786 }
787 
788 
789 
790 template <int dim, typename Number>
791 std::size_t
793 {
794  std::size_t total_memory_consumption = 0;
795  for (const auto &el : matrix_free_transfer_vector)
796  total_memory_consumption += el.memory_consumption();
797  return total_memory_consumption;
798 }
799 
800 
801 
802 // explicit instantiation
803 #include "mg_transfer_matrix_free.inst"
804 
805 
Utilities::pow
constexpr T pow(const T base, const int iexp)
Definition: utilities.h:476
LinearAlgebra::distributed::Vector< Number >
vectorization.h
MGLevelGlobalTransfer
Definition: mg_transfer.h:279
MGConstrainedDoFs
Definition: mg_constrained_dofs.h:45
StandardExceptions::ExcNotImplemented
static ::ExceptionBase & ExcNotImplemented()
MGTransferBlockMatrixFree::initialize_constraints
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
Definition: mg_transfer_matrix_free.cc:678
tria.h
tria_iterator.h
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
mg_transfer_internal.h
VectorizedArray
Definition: vectorization.h:395
MGTransferMatrixFree::do_restrict_add
void do_restrict_add(const unsigned int from_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
Definition: mg_transfer_matrix_free.cc:516
DoFHandler::n_components
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
VectorOperation::add
@ add
Definition: vector_operation.h:53
MGTransferBlockMatrixFree::restrict_and_add
virtual void restrict_and_add(const unsigned int from_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const override
Definition: mg_transfer_matrix_free.cc:768
MGTransferBlockMatrixFree::prolongate
virtual void prolongate(const unsigned int to_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const override
Definition: mg_transfer_matrix_free.cc:744
second
Point< 2 > second
Definition: grid_out.cc:4353
MGTransferBlockMatrixFree::memory_consumption
std::size_t memory_consumption() const
Definition: mg_transfer_matrix_free.cc:792
DoFHandler
Definition: dof_handler.h:205
LinearAlgebra::distributed::Vector::update_ghost_values
void update_ghost_values() const
la_parallel_vector.h
MatrixFreeOperators::BlockHelper::n_blocks
std::enable_if< IsBlockVector< VectorType >::value, unsigned int >::type n_blocks(const VectorType &vector)
Definition: operators.h:49
level
unsigned int level
Definition: grid_out.cc:4355
LinearAlgebra::distributed::Vector::get_partitioner
const std::shared_ptr< const Utilities::MPI::Partitioner > & get_partitioner() const
internal::MGTransfer::ElementInfo::element_is_continuous
bool element_is_continuous
Definition: mg_transfer_internal.h:91
mg_tools.h
StandardExceptions::ExcIndexRange
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
MGTransferMatrixFree::build
void build(const DoFHandler< dim, dim > &dof_handler, const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner >> &external_partitioners=std::vector< std::shared_ptr< const Utilities::MPI::Partitioner >>())
Definition: mg_transfer_matrix_free.cc:97
StandardExceptions::ExcMessage
static ::ExceptionBase & ExcMessage(std::string arg1)
TrilinosWrappers::internal::begin
VectorType::value_type * begin(VectorType &V)
Definition: trilinos_sparse_matrix.cc:51
MGTransferBlockMatrixFree::matrix_free_transfer_vector
std::vector< MGTransferMatrixFree< dim, Number > > matrix_free_transfer_vector
Definition: mg_transfer_matrix_free.h:480
DoFHandler::get_triangulation
const Triangulation< dim, spacedim > & get_triangulation() const
internal::MGTransfer::ElementInfo
Definition: mg_transfer_internal.h:79
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
BlockVectorBase< Vector< Number > >::n_blocks
unsigned int n_blocks() const
MemoryConsumption::memory_consumption
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
Definition: memory_consumption.h:268
Physics::Elasticity::Kinematics::b
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
VectorizedArrayBase< VectorizedArray< Number, width >, 1 >::size
static constexpr std::size_t size()
Definition: vectorization.h:261
LinearAlgebra::distributed::BlockVector
Definition: la_parallel_block_vector.h:85
MGTransferBlockMatrixFree::MGTransferBlockMatrixFree
MGTransferBlockMatrixFree()=default
fe.h
LinearAlgebra::distributed::Vector::local_size
size_type local_size() const
internal::MGTransfer::ElementInfo::n_components
unsigned int n_components
Definition: mg_transfer_internal.h:96
internal::evaluate_general
@ evaluate_general
Definition: tensor_product_kernels.h:42
AssertDimension
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
evaluation_kernels.h
LinearAlgebra::distributed::Vector::zero_out_ghosts
void zero_out_ghosts() const
MGTransferMatrixFree::restrict_and_add
virtual void restrict_and_add(const unsigned int from_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const override
Definition: mg_transfer_matrix_free.cc:264
MGTransferMatrixFree::MGTransferMatrixFree
MGTransferMatrixFree()
Definition: mg_transfer_matrix_free.cc:43
dof_tools.h
internal::FEEvaluationImplBasisChange
Definition: evaluation_kernels.h:668
function.h
MGTransferBlockMatrixFree::clear
void clear()
Definition: mg_transfer_matrix_free.cc:712
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
internal::MGTransfer::ElementInfo::prolongation_matrix_1d
std::vector< Number > prolongation_matrix_1d
Definition: mg_transfer_internal.h:116
dof_accessor.h
GridTools::shift
void shift(const Tensor< 1, spacedim > &shift_vector, Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:817
LinearAlgebra::distributed::Vector::reinit
void reinit(const size_type size, const bool omit_zeroing_entries=false)
internal::MGTransfer::ElementInfo::n_child_cell_dofs
unsigned int n_child_cell_dofs
Definition: mg_transfer_internal.h:103
mg_transfer_matrix_free.h
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
MGTransferMatrixFree::memory_consumption
std::size_t memory_consumption() const
Definition: mg_transfer_matrix_free.cc:639
BlockVectorBase< Vector< Number > >::block
BlockType & block(const unsigned int i)
logstream.h
first
Point< 2 > first
Definition: grid_out.cc:4352
LinearAlgebra::distributed::Vector::local_element
Number local_element(const size_type local_index) const
internal::MGTransfer::ElementInfo::fe_degree
unsigned int fe_degree
Definition: mg_transfer_internal.h:85
MGTransferMatrixFree::do_prolongate_add
void do_prolongate_add(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
Definition: mg_transfer_matrix_free.cc:377
MGTransferMatrixFree::initialize_constraints
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
Definition: mg_transfer_matrix_free.cc:67
MGTransferMatrixFree::prolongate
virtual void prolongate(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const override
Definition: mg_transfer_matrix_free.cc:183
MGTransferBlockMatrixFree::build
void build(const DoFHandler< dim, dim > &dof_handler)
Definition: mg_transfer_matrix_free.cc:721
MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > >::mg_constrained_dofs
SmartPointer< const MGConstrainedDoFs, MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > > > mg_constrained_dofs
Definition: mg_transfer.h:618
MGTransferMatrixFree::clear
void clear()
Definition: mg_transfer_matrix_free.cc:77