Reference documentation for deal.II version 9.0.0
mg_transfer_matrix_free.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2016 - 2018 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 #include <deal.II/base/vectorization.h>
20 #include <deal.II/lac/la_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/dofs/dof_accessor.h>
26 #include <deal.II/multigrid/mg_tools.h>
27 #include <deal.II/multigrid/mg_transfer_matrix_free.h>
28 #include <deal.II/multigrid/mg_transfer_internal.h>
29 
30 #include <deal.II/matrix_free/evaluation_kernels.h>
31 
32 #include <algorithm>
33 
34 DEAL_II_NAMESPACE_OPEN
35 
36 
37 template <int dim, typename Number>
39  :
40  fe_degree(0),
41  element_is_continuous(false),
42  n_components(0),
43  n_child_cell_dofs(0)
44 {}
45 
46 
47 
48 template <int dim, typename Number>
50  :
51  fe_degree(0),
52  element_is_continuous(false),
53  n_components(0),
54  n_child_cell_dofs(0)
55 {
56  this->mg_constrained_dofs = &mg_c;
57 }
58 
59 
60 
61 template <int dim, typename Number>
63 (const MGConstrainedDoFs &mg_c)
64 {
65  this->mg_constrained_dofs = &mg_c;
66 }
67 
68 
69 
70 template <int dim, typename Number>
72 {
74  fe_degree = 0;
75  element_is_continuous = false;
76  n_components = 0;
77  n_child_cell_dofs = 0;
78  level_dof_indices.clear();
79  parent_child_connect.clear();
80  dirichlet_indices.clear();
81  n_owned_level_cells.clear();
82  prolongation_matrix_1d.clear();
83  evaluation_data.clear();
84  weights_on_refined.clear();
85 }
86 
87 
88 template <int dim, typename Number>
90 (const DoFHandler<dim,dim> &mg_dof)
91 {
92  this->fill_and_communicate_copy_indices(mg_dof);
93 
94  std::vector<std::vector<Number> > weights_unvectorized;
95 
97 
98  internal::MGTransfer::setup_transfer<dim,Number>(mg_dof,
99  this->mg_constrained_dofs,
100  elem_info,
101  level_dof_indices,
102  parent_child_connect,
103  n_owned_level_cells,
104  dirichlet_indices,
105  weights_unvectorized,
106  this->copy_indices_global_mine,
107  this->ghosted_level_vector);
108  // unpack element info data
109  fe_degree = elem_info.fe_degree;
110  element_is_continuous = elem_info.element_is_continuous;
111  n_components = elem_info.n_components;
112  n_child_cell_dofs = elem_info.n_child_cell_dofs;
113 
114  // duplicate and put into vectorized array
115  prolongation_matrix_1d.resize(elem_info.prolongation_matrix_1d.size());
116  for (unsigned int i=0; i<elem_info.prolongation_matrix_1d.size(); i++)
117  prolongation_matrix_1d[i] = elem_info.prolongation_matrix_1d[i];
118 
119  // reshuffle into aligned vector of vectorized arrays
120  const unsigned int vec_size = VectorizedArray<Number>::n_array_elements;
121  const unsigned int n_levels = mg_dof.get_triangulation().n_global_levels();
122 
123  const unsigned int n_weights_per_cell = Utilities::fixed_power<dim>(3);
124  weights_on_refined.resize(n_levels-1);
125  for (unsigned int level = 1; level<n_levels; ++level)
126  {
127  weights_on_refined[level-1].resize(((n_owned_level_cells[level-1]+vec_size-1)/vec_size)*n_weights_per_cell);
128 
129  for (unsigned int c=0; c<n_owned_level_cells[level-1]; ++c)
130  {
131  const unsigned int comp = c/vec_size;
132  const unsigned int v = c%vec_size;
133  for (unsigned int i = 0; i<n_weights_per_cell; ++i)
134  {
135 
136  weights_on_refined[level-1][comp*n_weights_per_cell+i][v] = weights_unvectorized[level-1][c*n_weights_per_cell+i];
137  }
138  }
139  }
140 
141  evaluation_data.resize(n_child_cell_dofs);
142 }
143 
144 
145 
146 template <int dim, typename Number>
148 ::prolongate (const unsigned int to_level,
151 {
152  Assert ((to_level >= 1) && (to_level<=level_dof_indices.size()),
153  ExcIndexRange (to_level, 1, level_dof_indices.size()+1));
154 
155  AssertDimension(this->ghosted_level_vector[to_level].local_size(),
156  dst.local_size());
157  AssertDimension(this->ghosted_level_vector[to_level-1].local_size(),
158  src.local_size());
159 
160  this->ghosted_level_vector[to_level-1].copy_locally_owned_data_from(src);
161  this->ghosted_level_vector[to_level-1].update_ghost_values();
162  this->ghosted_level_vector[to_level] = 0.;
163 
164  // the implementation in do_prolongate_add is templated in the degree of the
165  // element (for efficiency reasons), so we need to find the appropriate
166  // kernel here...
167  if (fe_degree == 0)
168  do_prolongate_add<0>(to_level, this->ghosted_level_vector[to_level],
169  this->ghosted_level_vector[to_level-1]);
170  else if (fe_degree == 1)
171  do_prolongate_add<1>(to_level, this->ghosted_level_vector[to_level],
172  this->ghosted_level_vector[to_level-1]);
173  else if (fe_degree == 2)
174  do_prolongate_add<2>(to_level, this->ghosted_level_vector[to_level],
175  this->ghosted_level_vector[to_level-1]);
176  else if (fe_degree == 3)
177  do_prolongate_add<3>(to_level, this->ghosted_level_vector[to_level],
178  this->ghosted_level_vector[to_level-1]);
179  else if (fe_degree == 4)
180  do_prolongate_add<4>(to_level, this->ghosted_level_vector[to_level],
181  this->ghosted_level_vector[to_level-1]);
182  else if (fe_degree == 5)
183  do_prolongate_add<5>(to_level, this->ghosted_level_vector[to_level],
184  this->ghosted_level_vector[to_level-1]);
185  else if (fe_degree == 6)
186  do_prolongate_add<6>(to_level, this->ghosted_level_vector[to_level],
187  this->ghosted_level_vector[to_level-1]);
188  else if (fe_degree == 7)
189  do_prolongate_add<7>(to_level, this->ghosted_level_vector[to_level],
190  this->ghosted_level_vector[to_level-1]);
191  else if (fe_degree == 8)
192  do_prolongate_add<8>(to_level, this->ghosted_level_vector[to_level],
193  this->ghosted_level_vector[to_level-1]);
194  else if (fe_degree == 9)
195  do_prolongate_add<9>(to_level, this->ghosted_level_vector[to_level],
196  this->ghosted_level_vector[to_level-1]);
197  else if (fe_degree == 10)
198  do_prolongate_add<10>(to_level, this->ghosted_level_vector[to_level],
199  this->ghosted_level_vector[to_level-1]);
200  else
201  do_prolongate_add<-1>(to_level, this->ghosted_level_vector[to_level],
202  this->ghosted_level_vector[to_level-1]);
203 
204  this->ghosted_level_vector[to_level].compress(VectorOperation::add);
205  dst.copy_locally_owned_data_from(this->ghosted_level_vector[to_level]);
206 }
207 
208 
209 
210 template <int dim, typename Number>
212 ::restrict_and_add (const unsigned int from_level,
215 {
216  Assert ((from_level >= 1) && (from_level<=level_dof_indices.size()),
217  ExcIndexRange (from_level, 1, level_dof_indices.size()+1));
218 
219  AssertDimension(this->ghosted_level_vector[from_level].local_size(),
220  src.local_size());
221  AssertDimension(this->ghosted_level_vector[from_level-1].local_size(),
222  dst.local_size());
223 
224  this->ghosted_level_vector[from_level].copy_locally_owned_data_from(src);
225  this->ghosted_level_vector[from_level].update_ghost_values();
226  this->ghosted_level_vector[from_level-1] = 0.;
227 
228  if (fe_degree == 0)
229  do_restrict_add<0>(from_level, this->ghosted_level_vector[from_level-1],
230  this->ghosted_level_vector[from_level]);
231  else if (fe_degree == 1)
232  do_restrict_add<1>(from_level, this->ghosted_level_vector[from_level-1],
233  this->ghosted_level_vector[from_level]);
234  else if (fe_degree == 2)
235  do_restrict_add<2>(from_level, this->ghosted_level_vector[from_level-1],
236  this->ghosted_level_vector[from_level]);
237  else if (fe_degree == 3)
238  do_restrict_add<3>(from_level, this->ghosted_level_vector[from_level-1],
239  this->ghosted_level_vector[from_level]);
240  else if (fe_degree == 4)
241  do_restrict_add<4>(from_level, this->ghosted_level_vector[from_level-1],
242  this->ghosted_level_vector[from_level]);
243  else if (fe_degree == 5)
244  do_restrict_add<5>(from_level, this->ghosted_level_vector[from_level-1],
245  this->ghosted_level_vector[from_level]);
246  else if (fe_degree == 6)
247  do_restrict_add<6>(from_level, this->ghosted_level_vector[from_level-1],
248  this->ghosted_level_vector[from_level]);
249  else if (fe_degree == 7)
250  do_restrict_add<7>(from_level, this->ghosted_level_vector[from_level-1],
251  this->ghosted_level_vector[from_level]);
252  else if (fe_degree == 8)
253  do_restrict_add<8>(from_level, this->ghosted_level_vector[from_level-1],
254  this->ghosted_level_vector[from_level]);
255  else if (fe_degree == 9)
256  do_restrict_add<9>(from_level, this->ghosted_level_vector[from_level-1],
257  this->ghosted_level_vector[from_level]);
258  else if (fe_degree == 10)
259  do_restrict_add<10>(from_level, this->ghosted_level_vector[from_level-1],
260  this->ghosted_level_vector[from_level]);
261  else
262  // go to the non-templated version of the evaluator
263  do_restrict_add<-1>(from_level, this->ghosted_level_vector[from_level-1],
264  this->ghosted_level_vector[from_level]);
265 
266  this->ghosted_level_vector[from_level-1].compress(VectorOperation::add);
267  dst += this->ghosted_level_vector[from_level-1];
268 }
269 
270 
271 
272 namespace
273 {
274  template <int dim, int degree, typename Number>
275  void weight_dofs_on_child (const VectorizedArray<Number> *weights,
276  const unsigned int n_components,
277  const unsigned int fe_degree,
279  {
280  Assert(fe_degree > 0, ExcNotImplemented());
281  Assert(fe_degree < 100, ExcNotImplemented());
282  const int loop_length = degree != -1 ? 2*degree+1 : 2*fe_degree+1;
283  unsigned int degree_to_3 [100];
284  degree_to_3[0] = 0;
285  for (int i=1; i<loop_length-1; ++i)
286  degree_to_3[i] = 1;
287  degree_to_3[loop_length-1] = 2;
288  for (unsigned int c=0; c<n_components; ++c)
289  for (int k=0; k<(dim>2 ? loop_length : 1); ++k)
290  for (int j=0; j<(dim>1 ? loop_length : 1); ++j)
291  {
292  const unsigned int shift = 9*degree_to_3[k] + 3*degree_to_3[j];
293  data[0] *= weights[shift];
294  // loop bound as int avoids compiler warnings in case loop_length
295  // == 1 (polynomial degree 0)
296  for (int i=1; i<loop_length-1; ++i)
297  data[i] *= weights[shift+1];
298  data[loop_length-1] *= weights[shift+2];
299  data += loop_length;
300  }
301  }
302 }
303 
304 
305 
306 template <int dim, typename Number>
307 template <int degree>
309 ::do_prolongate_add (const unsigned int to_level,
312 {
313  const unsigned int vec_size = VectorizedArray<Number>::n_array_elements;
314  const unsigned int degree_size = (degree > -1 ? degree : fe_degree) + 1;
315  const unsigned int n_child_dofs_1d = 2*degree_size - element_is_continuous;
316  const unsigned int n_scalar_cell_dofs = Utilities::fixed_power<dim>(n_child_dofs_1d);
317  const unsigned int three_to_dim = Utilities::fixed_int_power<3,dim>::value;
318 
319  for (unsigned int cell=0; cell < n_owned_level_cells[to_level-1];
320  cell += vec_size)
321  {
322  const unsigned int n_lanes = cell+vec_size > n_owned_level_cells[to_level-1] ?
323  n_owned_level_cells[to_level-1] - cell : vec_size;
324 
325  // read from source vector
326  for (unsigned int v=0; v<n_lanes; ++v)
327  {
328  const unsigned int shift = internal::MGTransfer::compute_shift_within_children<dim>
329  (parent_child_connect[to_level-1][cell+v].second,
330  fe_degree+1-element_is_continuous, fe_degree);
331  const unsigned int *indices = &level_dof_indices[to_level-1][parent_child_connect[to_level-1][cell+v].first*n_child_cell_dofs+shift];
332  for (unsigned int c=0, m=0; c<n_components; ++c)
333  {
334  for (unsigned int k=0; k<(dim>2 ? degree_size : 1); ++k)
335  for (unsigned int j=0; j<(dim>1 ? degree_size : 1); ++j)
336  for (unsigned int i=0; i<degree_size; ++i, ++m)
337  evaluation_data[m][v] =
338  src.local_element(indices[c*n_scalar_cell_dofs +
339  k*n_child_dofs_1d*n_child_dofs_1d+
340  j*n_child_dofs_1d+i]);
341 
342  // apply Dirichlet boundary conditions on parent cell
343  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)
344  evaluation_data[*i][v] = 0.;
345  }
346  }
347 
348  AssertDimension(prolongation_matrix_1d.size(),
349  degree_size * n_child_dofs_1d);
350  // perform tensorized operation
351  if (element_is_continuous)
352  {
353  // must go through the components backwards because we want to write
354  // the output to the same array as the input
355  for (int c=n_components-1; c>=0; --c)
358  ::do_forward(prolongation_matrix_1d,
359  evaluation_data.begin()+ c*Utilities::fixed_power<dim>(degree_size),
360  evaluation_data.begin()+ c*n_scalar_cell_dofs,
361  fe_degree+1, 2*fe_degree+1);
362  weight_dofs_on_child<dim,degree,Number>(&weights_on_refined[to_level-1][(cell/vec_size)*three_to_dim],
363  n_components, fe_degree,
364  evaluation_data.begin());
365  }
366  else
367  {
368  for (int c=n_components-1; c>=0; --c)
371  ::do_forward(prolongation_matrix_1d,
372  evaluation_data.begin() + c*Utilities::fixed_power<dim>(degree_size),
373  evaluation_data.begin() + c*n_scalar_cell_dofs,
374  fe_degree+1, 2*fe_degree+2);
375  }
376 
377  // write into dst vector
378  const unsigned int *indices = &level_dof_indices[to_level][cell*
379  n_child_cell_dofs];
380  for (unsigned int v=0; v<n_lanes; ++v)
381  {
382  for (unsigned int i=0; i<n_child_cell_dofs; ++i)
383  dst.local_element(indices[i]) += evaluation_data[i][v];
384  indices += n_child_cell_dofs;
385  }
386  }
387 }
388 
389 
390 
391 template <int dim, typename Number>
392 template <int degree>
394 ::do_restrict_add (const unsigned int from_level,
397 {
398  const unsigned int vec_size = VectorizedArray<Number>::n_array_elements;
399  const unsigned int degree_size = (degree > -1 ? degree : fe_degree) + 1;
400  const unsigned int n_child_dofs_1d = 2*degree_size - element_is_continuous;
401  const unsigned int n_scalar_cell_dofs = Utilities::fixed_power<dim>(n_child_dofs_1d);
402  const unsigned int three_to_dim = Utilities::fixed_int_power<3,dim>::value;
403 
404  for (unsigned int cell=0; cell < n_owned_level_cells[from_level-1];
405  cell += vec_size)
406  {
407  const unsigned int n_lanes = cell+vec_size > n_owned_level_cells[from_level-1] ?
408  n_owned_level_cells[from_level-1] - cell : vec_size;
409 
410  // read from source vector
411  {
412  const unsigned int *indices = &level_dof_indices[from_level][cell*
413  n_child_cell_dofs];
414  for (unsigned int v=0; v<n_lanes; ++v)
415  {
416  for (unsigned int i=0; i<n_child_cell_dofs; ++i)
417  evaluation_data[i][v] = src.local_element(indices[i]);
418  indices += n_child_cell_dofs;
419  }
420  }
421 
422  AssertDimension(prolongation_matrix_1d.size(),
423  degree_size * n_child_dofs_1d);
424  // perform tensorized operation
425  if (element_is_continuous)
426  {
427  weight_dofs_on_child<dim,degree,Number>(&weights_on_refined[from_level-1][(cell/vec_size)*three_to_dim],
428  n_components, fe_degree,
429  &evaluation_data[0]);
430  for (unsigned int c=0; c<n_components; ++c)
433  ::do_backward(prolongation_matrix_1d, false,
434  evaluation_data.begin() + c*n_scalar_cell_dofs,
435  evaluation_data.begin() + c*Utilities::fixed_power<dim>(degree_size),
436  fe_degree+1, 2*fe_degree+1);
437  }
438  else
439  {
440  for (unsigned int c=0; c<n_components; ++c)
443  ::do_backward(prolongation_matrix_1d, false,
444  evaluation_data.begin() + c*n_scalar_cell_dofs,
445  evaluation_data.begin() + c*Utilities::fixed_power<dim>(degree_size),
446  fe_degree+1, 2*fe_degree+2);
447  }
448 
449  // write into dst vector
450  for (unsigned int v=0; v<n_lanes; ++v)
451  {
452  const unsigned int shift = internal::MGTransfer::compute_shift_within_children<dim>
453  (parent_child_connect[from_level-1][cell+v].second,
454  fe_degree+1-element_is_continuous, fe_degree);
455  AssertIndexRange(parent_child_connect[from_level-1][cell+v].first*
456  n_child_cell_dofs+n_child_cell_dofs-1,
457  level_dof_indices[from_level-1].size());
458  const unsigned int *indices = &level_dof_indices[from_level-1][parent_child_connect[from_level-1][cell+v].first*n_child_cell_dofs+shift];
459  for (unsigned int c=0, m=0; c<n_components; ++c)
460  {
461  // apply Dirichlet boundary conditions on parent cell
462  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)
463  evaluation_data[*i][v] = 0.;
464 
465  for (unsigned int k=0; k<(dim>2 ? degree_size : 1); ++k)
466  for (unsigned int j=0; j<(dim>1 ? degree_size : 1); ++j)
467  for (unsigned int i=0; i<degree_size; ++i, ++m)
468  dst.local_element(indices[c*n_scalar_cell_dofs +
469  k*n_child_dofs_1d*n_child_dofs_1d+
470  j*n_child_dofs_1d+i])
471  += evaluation_data[m][v];
472  }
473  }
474  }
475 }
476 
477 
478 
479 template <int dim, typename Number>
480 std::size_t
482 {
483  std::size_t memory = MGLevelGlobalTransfer<LinearAlgebra::distributed::Vector<Number> >::memory_consumption();
484  memory += MemoryConsumption::memory_consumption(level_dof_indices);
485  memory += MemoryConsumption::memory_consumption(parent_child_connect);
486  memory += MemoryConsumption::memory_consumption(n_owned_level_cells);
487  memory += MemoryConsumption::memory_consumption(prolongation_matrix_1d);
488  memory += MemoryConsumption::memory_consumption(evaluation_data);
489  memory += MemoryConsumption::memory_consumption(weights_on_refined);
490  memory += MemoryConsumption::memory_consumption(dirichlet_indices);
491  return memory;
492 }
493 
494 
495 
496 template <int dim, typename Number>
498  : same_for_all (true)
499 {
500  matrix_free_transfer_vector.emplace_back(mg_c);
501 }
502 
503 
504 
505 template <int dim, typename Number>
507 (const std::vector<MGConstrainedDoFs> &mg_c)
508  : same_for_all (false)
509 {
510  for (unsigned int i=0; i<mg_c.size(); ++i)
511  matrix_free_transfer_vector.emplace_back(mg_c[i]);
512 }
513 
514 
515 
516 template <int dim, typename Number>
518 (const MGConstrainedDoFs &mg_c)
519 {
520  Assert(same_for_all,
521  ExcMessage("This object was initialized with support for usage with "
522  "one DoFHandler for each block, but this method assumes "
523  "that the same DoFHandler is used for all the blocks!"));
524  AssertDimension(matrix_free_transfer_vector.size(),1);
525 
526  matrix_free_transfer_vector[0].initialize_constraints(mg_c);
527 }
528 
529 
530 
531 template <int dim, typename Number>
533 (const std::vector<MGConstrainedDoFs> &mg_c)
534 {
535  Assert(!same_for_all,
536  ExcMessage("This object was initialized with support for using "
537  "the same DoFHandler for all the blocks, but this "
538  "method assumes that there is a separate DoFHandler "
539  "for each block!"));
540  AssertDimension(matrix_free_transfer_vector.size(),mg_c.size());
541 
542  for (unsigned int i = 0; i<mg_c.size(); ++i)
543  matrix_free_transfer_vector[i].initialize_constraints(mg_c[i]);
544 }
545 
546 
547 
548 template <int dim, typename Number>
550 {
551  matrix_free_transfer_vector.clear();
552 }
553 
554 
555 
556 template <int dim, typename Number>
558 (const DoFHandler<dim,dim> &mg_dof)
559 {
560  AssertDimension(matrix_free_transfer_vector.size(),1);
561  matrix_free_transfer_vector[0].build(mg_dof);
562 }
563 
564 
565 
566 template <int dim, typename Number>
568 (const std::vector<const DoFHandler<dim,dim>*> &mg_dof)
569 {
570  AssertDimension (matrix_free_transfer_vector.size(),mg_dof.size());
571  for (unsigned int i=0; i<mg_dof.size(); ++i)
572  matrix_free_transfer_vector[i].build(*mg_dof[i]);
573 }
574 
575 
576 
577 template <int dim, typename Number>
579 ::prolongate (const unsigned int to_level,
582 {
583  const unsigned int n_blocks = src.n_blocks();
584  AssertDimension(dst.n_blocks(), n_blocks);
585 
586  if (!same_for_all)
587  AssertDimension (matrix_free_transfer_vector.size(), n_blocks);
588 
589  for (unsigned int b = 0; b < n_blocks; ++b)
590  {
591  const unsigned int data_block = same_for_all ? 0 : b;
592  matrix_free_transfer_vector[data_block].prolongate(to_level, dst.block(b), src.block(b));
593  }
594 }
595 
596 
597 
598 template <int dim, typename Number>
600 ::restrict_and_add (const unsigned int from_level,
603 {
604  const unsigned int n_blocks = src.n_blocks();
605  AssertDimension(dst.n_blocks(), n_blocks);
606 
607  if (!same_for_all)
608  AssertDimension (matrix_free_transfer_vector.size(), n_blocks);
609 
610  for (unsigned int b = 0; b < n_blocks; ++b)
611  {
612  const unsigned int data_block = same_for_all ? 0 : b;
613  matrix_free_transfer_vector[data_block].restrict_and_add(from_level, dst.block(b), src.block(b));
614  }
615 }
616 
617 
618 
619 template <int dim, typename Number>
620 std::size_t
622 {
623  std::size_t total_memory_consumption= 0;
624  for (const auto &el: matrix_free_transfer_vector)
625  total_memory_consumption += el.memory_consumption();
626  return total_memory_consumption;
627 }
628 
629 
630 
631 // explicit instantiation
632 #include "mg_transfer_matrix_free.inst"
633 
634 
635 DEAL_II_NAMESPACE_CLOSE
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:1248
virtual void restrict_and_add(const unsigned int from_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
void copy_locally_owned_data_from(const Vector< Number2 > &src)
std::size_t memory_consumption() const
void build(const DoFHandler< dim, dim > &mg_dof)
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
#define AssertIndexRange(index, range)
Definition: exceptions.h:1284
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
std::size_t memory_consumption() const
std::vector< MGTransferMatrixFree< dim, Number > > matrix_free_transfer_vector
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1142
SmartPointer< const MGConstrainedDoFs, MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > > > mg_constrained_dofs
Definition: mg_transfer.h:486
virtual void prolongate(const unsigned int to_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const
void do_prolongate_add(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
MGTransferBlockMatrixFree()=default
Number local_element(const size_type local_index) const
virtual void restrict_and_add(const unsigned int from_level, LinearAlgebra::distributed::BlockVector< Number > &dst, const LinearAlgebra::distributed::BlockVector< Number > &src) const
static ::ExceptionBase & ExcNotImplemented()
virtual void prolongate(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
BlockType & block(const unsigned int i)
void build(const DoFHandler< dim, dim > &mg_dof)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)