Reference documentation for deal.II version Git a285ab78cc 2020-08-04 19:09:24 +0200
\(\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)
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;
85  level_dof_indices.clear();
86  parent_child_connect.clear();
87  dirichlet_indices.clear();
88  n_owned_level_cells.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 
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)
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,
125  weights_unvectorized,
128 
129  // reinit the ghosted level vector in case its partitioner differs from the
130  // one the user intends to use externally
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())
137  else
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 
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  }
214 
216  src_inplace ? src : this->ghosted_level_vector[to_level - 1];
218  dst_inplace ? dst : this->ghosted_level_vector[to_level];
219 
220  dst_vec = 0.;
221 
222  src_vec.update_ghost_values();
223  // the implementation in do_prolongate_add is templated in the degree of the
224  // element (for efficiency reasons), so we need to find the appropriate
225  // kernel here...
226  if (fe_degree == 0)
227  do_prolongate_add<0>(to_level, dst_vec, src_vec);
228  else if (fe_degree == 1)
229  do_prolongate_add<1>(to_level, dst_vec, src_vec);
230  else if (fe_degree == 2)
231  do_prolongate_add<2>(to_level, dst_vec, src_vec);
232  else if (fe_degree == 3)
233  do_prolongate_add<3>(to_level, dst_vec, src_vec);
234  else if (fe_degree == 4)
235  do_prolongate_add<4>(to_level, dst_vec, src_vec);
236  else if (fe_degree == 5)
237  do_prolongate_add<5>(to_level, dst_vec, src_vec);
238  else if (fe_degree == 6)
239  do_prolongate_add<6>(to_level, dst_vec, src_vec);
240  else if (fe_degree == 7)
241  do_prolongate_add<7>(to_level, dst_vec, src_vec);
242  else if (fe_degree == 8)
243  do_prolongate_add<8>(to_level, dst_vec, src_vec);
244  else if (fe_degree == 9)
245  do_prolongate_add<9>(to_level, dst_vec, src_vec);
246  else if (fe_degree == 10)
247  do_prolongate_add<10>(to_level, dst_vec, src_vec);
248  else
249  do_prolongate_add<-1>(to_level, dst_vec, src_vec);
250 
251  dst_vec.compress(VectorOperation::add);
252  if (dst_inplace == false)
254 
255  if (src_inplace == true)
256  src.zero_out_ghosts();
257 }
258 
259 
260 
261 template <int dim, typename Number>
262 void
264  const unsigned int from_level,
267 {
268  Assert((from_level >= 1) && (from_level <= level_dof_indices.size()),
269  ExcIndexRange(from_level, 1, level_dof_indices.size() + 1));
270 
271  const bool src_inplace =
272  src.get_partitioner().get() == this->vector_partitioners[from_level].get();
273  if (src_inplace == false)
274  {
275  if (this->ghosted_level_vector[from_level].get_partitioner().get() !=
276  this->vector_partitioners[from_level].get())
277  this->ghosted_level_vector[from_level].reinit(
278  this->vector_partitioners[from_level]);
279  this->ghosted_level_vector[from_level].copy_locally_owned_data_from(src);
280  }
281 
282  const bool dst_inplace = dst.get_partitioner().get() ==
283  this->vector_partitioners[from_level - 1].get();
284  if (dst_inplace == false)
285  {
286  if (this->ghosted_level_vector[from_level - 1].get_partitioner().get() !=
287  this->vector_partitioners[from_level - 1].get())
288  this->ghosted_level_vector[from_level - 1].reinit(
289  this->vector_partitioners[from_level - 1]);
290  AssertDimension(this->ghosted_level_vector[from_level - 1].local_size(),
291  dst.local_size());
292  this->ghosted_level_vector[from_level - 1] = 0.;
293  }
294 
296  src_inplace ? src : this->ghosted_level_vector[from_level];
298  dst_inplace ? dst : this->ghosted_level_vector[from_level - 1];
299 
300  src_vec.update_ghost_values();
301 
302  if (fe_degree == 0)
303  do_restrict_add<0>(from_level, dst_vec, src_vec);
304  else if (fe_degree == 1)
305  do_restrict_add<1>(from_level, dst_vec, src_vec);
306  else if (fe_degree == 2)
307  do_restrict_add<2>(from_level, dst_vec, src_vec);
308  else if (fe_degree == 3)
309  do_restrict_add<3>(from_level, dst_vec, src_vec);
310  else if (fe_degree == 4)
311  do_restrict_add<4>(from_level, dst_vec, src_vec);
312  else if (fe_degree == 5)
313  do_restrict_add<5>(from_level, dst_vec, src_vec);
314  else if (fe_degree == 6)
315  do_restrict_add<6>(from_level, dst_vec, src_vec);
316  else if (fe_degree == 7)
317  do_restrict_add<7>(from_level, dst_vec, src_vec);
318  else if (fe_degree == 8)
319  do_restrict_add<8>(from_level, dst_vec, src_vec);
320  else if (fe_degree == 9)
321  do_restrict_add<9>(from_level, dst_vec, src_vec);
322  else if (fe_degree == 10)
323  do_restrict_add<10>(from_level, dst_vec, src_vec);
324  else
325  // go to the non-templated version of the evaluator
326  do_restrict_add<-1>(from_level, dst_vec, src_vec);
327 
328  dst_vec.compress(VectorOperation::add);
329  if (dst_inplace == false)
330  dst += dst_vec;
331 
332  if (src_inplace == true)
333  src.zero_out_ghosts();
334 }
335 
336 
337 
338 namespace
339 {
340  template <int dim, int degree, typename Number>
341  void
342  weight_dofs_on_child(const VectorizedArray<Number> *weights,
343  const unsigned int n_components,
344  const unsigned int fe_degree,
346  {
347  Assert(fe_degree > 0, ExcNotImplemented());
348  Assert(fe_degree < 100, ExcNotImplemented());
349  const int loop_length = degree != -1 ? 2 * degree + 1 : 2 * fe_degree + 1;
350  unsigned int degree_to_3[100];
351  degree_to_3[0] = 0;
352  for (int i = 1; i < loop_length - 1; ++i)
353  degree_to_3[i] = 1;
354  degree_to_3[loop_length - 1] = 2;
355  for (unsigned int c = 0; c < n_components; ++c)
356  for (int k = 0; k < (dim > 2 ? loop_length : 1); ++k)
357  for (int j = 0; j < (dim > 1 ? loop_length : 1); ++j)
358  {
359  const unsigned int shift = 9 * degree_to_3[k] + 3 * degree_to_3[j];
360  data[0] *= weights[shift];
361  // loop bound as int avoids compiler warnings in case loop_length
362  // == 1 (polynomial degree 0)
363  for (int i = 1; i < loop_length - 1; ++i)
364  data[i] *= weights[shift + 1];
365  data[loop_length - 1] *= weights[shift + 2];
366  data += loop_length;
367  }
368  }
369 } // namespace
370 
371 
372 
373 template <int dim, typename Number>
374 template <int degree>
375 void
377  const unsigned int to_level,
380 {
381  const unsigned int vec_size = VectorizedArray<Number>::size();
382  const unsigned int degree_size = (degree > -1 ? degree : fe_degree) + 1;
383  const unsigned int n_child_dofs_1d = 2 * degree_size - element_is_continuous;
384  const unsigned int n_scalar_cell_dofs =
385  Utilities::fixed_power<dim>(n_child_dofs_1d);
386  constexpr unsigned int three_to_dim = Utilities::pow(3, dim);
387 
388  // If we have user defined MG constraints, we must create
389  // a non-const, ghosted version of the source vector to distribute
390  // constraints.
391  const LinearAlgebra::distributed::Vector<Number> *to_use = &src;
393  if (this->mg_constrained_dofs != nullptr &&
395  .get_local_lines()
396  .size() > 0)
397  {
399 
400  // Distribute any user defined constraints
402  .distribute(copy_src);
403 
404  // Re-initialize new ghosted vector with correct constraints
405  new_src.reinit(copy_src);
406  new_src = copy_src;
407  new_src.update_ghost_values();
408 
409  to_use = &new_src;
410  }
411 
412  for (unsigned int cell = 0; cell < n_owned_level_cells[to_level - 1];
413  cell += vec_size)
414  {
415  const unsigned int n_lanes =
416  cell + vec_size > n_owned_level_cells[to_level - 1] ?
417  n_owned_level_cells[to_level - 1] - cell :
418  vec_size;
419 
420  // read from source vector
421  for (unsigned int v = 0; v < n_lanes; ++v)
422  {
423  const unsigned int shift =
424  internal::MGTransfer::compute_shift_within_children<dim>(
425  parent_child_connect[to_level - 1][cell + v].second,
427  fe_degree);
428  const unsigned int *indices =
429  &level_dof_indices[to_level - 1]
430  [parent_child_connect[to_level - 1][cell + v]
431  .first *
433  shift];
434  for (unsigned int c = 0, m = 0; c < n_components; ++c)
435  {
436  for (unsigned int k = 0; k < (dim > 2 ? degree_size : 1); ++k)
437  for (unsigned int j = 0; j < (dim > 1 ? degree_size : 1); ++j)
438  for (unsigned int i = 0; i < degree_size; ++i, ++m)
439  evaluation_data[m][v] = to_use->local_element(
440  indices[c * n_scalar_cell_dofs +
441  k * n_child_dofs_1d * n_child_dofs_1d +
442  j * n_child_dofs_1d + i]);
443 
444  // apply Dirichlet boundary conditions on parent cell
445  for (std::vector<unsigned short>::const_iterator i =
446  dirichlet_indices[to_level - 1][cell + v].begin();
447  i != dirichlet_indices[to_level - 1][cell + v].end();
448  ++i)
449  evaluation_data[*i][v] = 0.;
450  }
451  }
452 
454  degree_size * n_child_dofs_1d);
455  // perform tensorized operation
456  if (element_is_continuous)
457  {
458  // must go through the components backwards because we want to write
459  // the output to the same array as the input
460  for (int c = n_components - 1; c >= 0; --c)
462  dim,
463  degree + 1,
464  2 * degree + 1,
465  1,
468  do_forward(prolongation_matrix_1d,
470  c * Utilities::fixed_power<dim>(degree_size),
471  evaluation_data.begin() + c * n_scalar_cell_dofs,
472  fe_degree + 1,
473  2 * fe_degree + 1);
474  weight_dofs_on_child<dim, degree, Number>(
475  &weights_on_refined[to_level - 1][(cell / vec_size) * three_to_dim],
476  n_components,
477  fe_degree,
479  }
480  else
481  {
482  for (int c = n_components - 1; c >= 0; --c)
484  dim,
485  degree + 1,
486  2 * degree + 2,
487  1,
490  do_forward(prolongation_matrix_1d,
492  c * Utilities::fixed_power<dim>(degree_size),
493  evaluation_data.begin() + c * n_scalar_cell_dofs,
494  fe_degree + 1,
495  2 * fe_degree + 2);
496  }
497 
498  // write into dst vector
499  const unsigned int *indices =
500  &level_dof_indices[to_level][cell * n_child_cell_dofs];
501  for (unsigned int v = 0; v < n_lanes; ++v)
502  {
503  for (unsigned int i = 0; i < n_child_cell_dofs; ++i)
504  dst.local_element(indices[i]) += evaluation_data[i][v];
505  indices += n_child_cell_dofs;
506  }
507  }
508 }
509 
510 
511 
512 template <int dim, typename Number>
513 template <int degree>
514 void
516  const unsigned int from_level,
519 {
520  const unsigned int vec_size = VectorizedArray<Number>::size();
521  const unsigned int degree_size = (degree > -1 ? degree : fe_degree) + 1;
522  const unsigned int n_child_dofs_1d = 2 * degree_size - element_is_continuous;
523  const unsigned int n_scalar_cell_dofs =
524  Utilities::fixed_power<dim>(n_child_dofs_1d);
525  constexpr unsigned int three_to_dim = Utilities::pow(3, dim);
526 
527  for (unsigned int cell = 0; cell < n_owned_level_cells[from_level - 1];
528  cell += vec_size)
529  {
530  const unsigned int n_lanes =
531  cell + vec_size > n_owned_level_cells[from_level - 1] ?
532  n_owned_level_cells[from_level - 1] - cell :
533  vec_size;
534 
535  // read from source vector
536  {
537  const unsigned int *indices =
538  &level_dof_indices[from_level][cell * n_child_cell_dofs];
539  for (unsigned int v = 0; v < n_lanes; ++v)
540  {
541  for (unsigned int i = 0; i < n_child_cell_dofs; ++i)
542  evaluation_data[i][v] = src.local_element(indices[i]);
543  indices += n_child_cell_dofs;
544  }
545  }
546 
548  degree_size * n_child_dofs_1d);
549  // perform tensorized operation
550  if (element_is_continuous)
551  {
552  weight_dofs_on_child<dim, degree, Number>(
553  &weights_on_refined[from_level - 1]
554  [(cell / vec_size) * three_to_dim],
555  n_components,
556  fe_degree,
558  for (unsigned int c = 0; c < n_components; ++c)
560  dim,
561  degree + 1,
562  2 * degree + 1,
563  1,
566  do_backward(prolongation_matrix_1d,
567  false,
568  evaluation_data.begin() + c * n_scalar_cell_dofs,
570  c * Utilities::fixed_power<dim>(degree_size),
571  fe_degree + 1,
572  2 * fe_degree + 1);
573  }
574  else
575  {
576  for (unsigned int c = 0; c < n_components; ++c)
578  dim,
579  degree + 1,
580  2 * degree + 2,
581  1,
584  do_backward(prolongation_matrix_1d,
585  false,
586  evaluation_data.begin() + c * n_scalar_cell_dofs,
588  c * Utilities::fixed_power<dim>(degree_size),
589  fe_degree + 1,
590  2 * fe_degree + 2);
591  }
592 
593  // write into dst vector
594  for (unsigned int v = 0; v < n_lanes; ++v)
595  {
596  const unsigned int shift =
597  internal::MGTransfer::compute_shift_within_children<dim>(
598  parent_child_connect[from_level - 1][cell + v].second,
600  fe_degree);
602  parent_child_connect[from_level - 1][cell + v].first *
604  n_child_cell_dofs - 1,
605  level_dof_indices[from_level - 1].size());
606  const unsigned int *indices =
607  &level_dof_indices[from_level - 1]
608  [parent_child_connect[from_level - 1][cell + v]
609  .first *
611  shift];
612  for (unsigned int c = 0, m = 0; c < n_components; ++c)
613  {
614  // apply Dirichlet boundary conditions on parent cell
615  for (std::vector<unsigned short>::const_iterator i =
616  dirichlet_indices[from_level - 1][cell + v].begin();
617  i != dirichlet_indices[from_level - 1][cell + v].end();
618  ++i)
619  evaluation_data[*i][v] = 0.;
620 
621  for (unsigned int k = 0; k < (dim > 2 ? degree_size : 1); ++k)
622  for (unsigned int j = 0; j < (dim > 1 ? degree_size : 1); ++j)
623  for (unsigned int i = 0; i < degree_size; ++i, ++m)
624  dst.local_element(
625  indices[c * n_scalar_cell_dofs +
626  k * n_child_dofs_1d * n_child_dofs_1d +
627  j * n_child_dofs_1d + i]) +=
628  evaluation_data[m][v];
629  }
630  }
631  }
632 }
633 
634 
635 
636 template <int dim, typename Number>
637 std::size_t
639 {
640  std::size_t memory = MGLevelGlobalTransfer<
649  return memory;
650 }
651 
652 
653 
654 template <int dim, typename Number>
656  const MGConstrainedDoFs &mg_c)
657  : same_for_all(true)
658 {
659  matrix_free_transfer_vector.emplace_back(mg_c);
660 }
661 
662 
663 
664 template <int dim, typename Number>
666  const std::vector<MGConstrainedDoFs> &mg_c)
667  : same_for_all(false)
668 {
669  for (const auto &constrained_block_dofs : mg_c)
670  matrix_free_transfer_vector.emplace_back(constrained_block_dofs);
671 }
672 
673 
674 
675 template <int dim, typename Number>
676 void
678  const MGConstrainedDoFs &mg_c)
679 {
681  ExcMessage("This object was initialized with support for usage with "
682  "one DoFHandler for each block, but this method assumes "
683  "that the same DoFHandler is used for all the blocks!"));
685 
686  matrix_free_transfer_vector[0].initialize_constraints(mg_c);
687 }
688 
689 
690 
691 template <int dim, typename Number>
692 void
694  const std::vector<MGConstrainedDoFs> &mg_c)
695 {
697  ExcMessage("This object was initialized with support for using "
698  "the same DoFHandler for all the blocks, but this "
699  "method assumes that there is a separate DoFHandler "
700  "for each block!"));
701  AssertDimension(matrix_free_transfer_vector.size(), mg_c.size());
702 
703  for (unsigned int i = 0; i < mg_c.size(); ++i)
705 }
706 
707 
708 
709 template <int dim, typename Number>
710 void
712 {
714 }
715 
716 
717 
718 template <int dim, typename Number>
719 void
721  const DoFHandler<dim, dim> &dof_handler)
722 {
724  matrix_free_transfer_vector[0].build(dof_handler);
725 }
726 
727 
728 
729 template <int dim, typename Number>
730 void
732  const std::vector<const DoFHandler<dim, dim> *> &dof_handler)
733 {
734  AssertDimension(matrix_free_transfer_vector.size(), dof_handler.size());
735  for (unsigned int i = 0; i < dof_handler.size(); ++i)
736  matrix_free_transfer_vector[i].build(*dof_handler[i]);
737 }
738 
739 
740 
741 template <int dim, typename Number>
742 void
744  const unsigned int to_level,
747 {
748  const unsigned int n_blocks = src.n_blocks();
750 
751  if (!same_for_all)
753 
754  for (unsigned int b = 0; b < n_blocks; ++b)
755  {
756  const unsigned int data_block = same_for_all ? 0 : b;
757  matrix_free_transfer_vector[data_block].prolongate(to_level,
758  dst.block(b),
759  src.block(b));
760  }
761 }
762 
763 
764 
765 template <int dim, typename Number>
766 void
768  const unsigned int from_level,
771 {
772  const unsigned int n_blocks = src.n_blocks();
774 
775  if (!same_for_all)
777 
778  for (unsigned int b = 0; b < n_blocks; ++b)
779  {
780  const unsigned int data_block = same_for_all ? 0 : b;
781  matrix_free_transfer_vector[data_block].restrict_and_add(from_level,
782  dst.block(b),
783  src.block(b));
784  }
785 }
786 
787 
788 
789 template <int dim, typename Number>
790 std::size_t
792 {
793  std::size_t total_memory_consumption = 0;
794  for (const auto &el : matrix_free_transfer_vector)
795  total_memory_consumption += el.memory_consumption();
796  return total_memory_consumption;
797 }
798 
799 
800 
801 // explicit instantiation
802 #include "mg_transfer_matrix_free.inst"
803 
804 
unsigned int max_level() const
const Triangulation< dim, spacedim > & get_triangulation() const
void copy_locally_owned_data_from(const Vector< Number2, MemorySpace > &src)
std::vector< Table< 2, unsigned int > > copy_indices_global_mine
Definition: mg_transfer.h:558
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
Definition: matrix_block.h:618
void do_restrict_add(const unsigned int from_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1568
std::vector< unsigned int > n_owned_level_cells
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 >>())
AlignedVector< VectorizedArray< Number > > evaluation_data
std::size_t memory_consumption() const
void build(const DoFHandler< dim, dim > &dof_handler)
std::vector< std::vector< std::vector< unsigned short > > > dirichlet_indices
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
#define AssertIndexRange(index, range)
Definition: exceptions.h:1636
pointer data()
void reinit(const size_type size, const bool omit_zeroing_entries=false)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
void resize(const size_type size_in)
std::size_t memory_consumption() const
constexpr T pow(const T base, const int iexp)
Definition: utilities.h:461
size_type size() const
Definition: index_set.h:1632
Number local_element(const size_type local_index) const
static ::ExceptionBase & ExcMessage(std::string arg1)
AlignedVector< VectorizedArray< Number > > prolongation_matrix_1d
#define Assert(cond, exc)
Definition: exceptions.h:1411
virtual void restrict_and_add(const unsigned int from_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const override
virtual void restrict_and_add(const unsigned int from_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const override
virtual void prolongate(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const override
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:363
unsigned int level
Definition: grid_out.cc:4341
std::vector< AlignedVector< VectorizedArray< Number > > > weights_on_refined
std::enable_if< IsBlockVector< VectorType >::value, unsigned int >::type n_blocks(const VectorType &vector)
Definition: operators.h:49
void fill_and_communicate_copy_indices(const DoFHandler< dim, spacedim > &dof_handler)
std::vector< std::vector< unsigned int > > level_dof_indices
const IndexSet & get_local_lines() const
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
Point< 2 > first
Definition: grid_out.cc:4338
void do_prolongate_add(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
SmartPointer< const MGConstrainedDoFs, MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > > > mg_constrained_dofs
Definition: mg_transfer.h:607
void distribute(VectorType &vec) const
std::vector< std::vector< std::pair< unsigned int, unsigned int > > > parent_child_connect
virtual void prolongate(const unsigned int to_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const override
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:362
VectorType::value_type * begin(VectorType &V)
std::vector< MGTransferMatrixFree< dim, Number > > matrix_free_transfer_vector
const std::shared_ptr< const Utilities::MPI::Partitioner > & get_partitioner() const
MGTransferBlockMatrixFree()=default
iterator begin()
size_type size() const
static ::ExceptionBase & ExcNotImplemented()
void resize(const unsigned int new_minlevel, const unsigned int new_maxlevel)
MGLevelObject< std::shared_ptr< const Utilities::MPI::Partitioner > > vector_partitioners
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
const AffineConstraints< double > & get_user_constraint_matrix(const unsigned int level) const
BlockType & block(const unsigned int i)
MGLevelObject< LinearAlgebra::distributed::Vector< Number > > ghosted_level_vector
Definition: mg_transfer.h:627
void shift(const Tensor< 1, spacedim > &shift_vector, Triangulation< dim, spacedim > &triangulation)
Definition: grid_tools.cc:815
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)