deal.II version GIT relicensing-2167-g9622207b8f 2024-11-21 12:40:00+00:00
\(\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\}}\)
Loading...
Searching...
No Matches
mg_transfer_matrix_free.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2016 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15
19
22
23#include <deal.II/fe/fe.h>
24
25#include <deal.II/grid/tria.h>
27
29
32
35#include <deal.II/multigrid/mg_transfer_matrix_free.templates.h>
36
37#include <algorithm>
38
40
41
42template <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
52template <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
65template <int dim, typename Number>
66void
68 const MGConstrainedDoFs &mg_c)
69{
70 this->mg_constrained_dofs = &mg_c;
71}
72
73
74
75template <int dim, typename Number>
76void
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
96template <int dim, typename Number>
97void
99 const DoFHandler<dim> &dof_handler,
100 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
101 &external_partitioners)
102{
103 Assert(dof_handler.has_level_dofs(),
105 "The underlying DoFHandler object has not had its "
106 "distribute_mg_dofs() function called, but this is a prerequisite "
107 "for multigrid transfers. You will need to call this function, "
108 "probably close to where you already call distribute_dofs()."));
109 if (external_partitioners.size() > 0)
110 {
111 Assert(this->initialize_dof_vector == nullptr,
112 ExcMessage("An initialize_dof_vector function has already "
113 "been registered in the constructor!"));
114
115 this->initialize_dof_vector =
116 [external_partitioners](
117 const unsigned int level,
119 vec.reinit(external_partitioners[level]);
120 };
121 }
122
123 this->fill_and_communicate_copy_indices(dof_handler);
124
125 vector_partitioners.resize(0,
126 dof_handler.get_triangulation().n_global_levels() -
127 1);
128 for (unsigned int level = 0; level <= this->ghosted_level_vector.max_level();
129 ++level)
130 vector_partitioners[level] =
131 this->ghosted_level_vector[level].get_partitioner();
132
133 std::vector<std::vector<Number>> weights_unvectorized;
134
136
137 internal::MGTransfer::setup_transfer<dim, Number>(
138 dof_handler,
139 this->mg_constrained_dofs,
140 external_partitioners,
141 elem_info,
142 level_dof_indices,
143 parent_child_connect,
144 n_owned_level_cells,
145 dirichlet_indices,
146 weights_unvectorized,
147 this->copy_indices_global_mine,
148 vector_partitioners);
149
150 // reinit the ghosted level vector in case its partitioner differs from the
151 // one the user intends to use externally
152 this->ghosted_level_vector.resize(0, vector_partitioners.max_level());
153 for (unsigned int level = 0; level <= vector_partitioners.max_level();
154 ++level)
155 if (external_partitioners.size() == vector_partitioners.max_level() + 1 &&
156 external_partitioners[level].get() == vector_partitioners[level].get())
157 this->ghosted_level_vector[level].reinit(0);
158 else
159 this->ghosted_level_vector[level].reinit(vector_partitioners[level]);
160
161 // unpack element info data
162 fe_degree = elem_info.fe_degree;
163 element_is_continuous = elem_info.element_is_continuous;
164 n_components = elem_info.n_components;
165 n_child_cell_dofs = elem_info.n_child_cell_dofs;
166
167 // duplicate and put into vectorized array
168 prolongation_matrix_1d.resize(elem_info.prolongation_matrix_1d.size());
169 for (unsigned int i = 0; i < elem_info.prolongation_matrix_1d.size(); ++i)
170 prolongation_matrix_1d[i] = elem_info.prolongation_matrix_1d[i];
171
172 // reshuffle into aligned vector of vectorized arrays
173 const unsigned int vec_size = VectorizedArray<Number>::size();
174 const unsigned int n_levels =
175 dof_handler.get_triangulation().n_global_levels();
176
177 const unsigned int n_weights_per_cell = Utilities::fixed_power<dim>(3);
178 weights_on_refined.resize(n_levels - 1);
179 for (unsigned int level = 1; level < n_levels; ++level)
180 {
181 weights_on_refined[level - 1].resize(
182 ((n_owned_level_cells[level - 1] + vec_size - 1) / vec_size) *
183 n_weights_per_cell);
184
185 for (unsigned int c = 0; c < n_owned_level_cells[level - 1]; ++c)
186 {
187 const unsigned int comp = c / vec_size;
188 const unsigned int v = c % vec_size;
189 for (unsigned int i = 0; i < n_weights_per_cell; ++i)
190 {
191 weights_on_refined[level - 1][comp * n_weights_per_cell + i][v] =
192 weights_unvectorized[level - 1][c * n_weights_per_cell + i];
193 }
194 }
195 }
196
197 evaluation_data.resize(n_child_cell_dofs);
198}
199
200
201
202template <int dim, typename Number>
203void
205 const DoFHandler<dim> &dof_handler,
206 const std::function<void(const unsigned int,
208 &initialize_dof_vector)
209{
210 if (initialize_dof_vector)
211 {
212 const unsigned int n_levels =
213 dof_handler.get_triangulation().n_global_levels();
214
215 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
216 external_partitioners(n_levels);
217
218 for (unsigned int level = 0; level < n_levels; ++level)
219 {
221 initialize_dof_vector(level, vector);
222 external_partitioners[level] = vector.get_partitioner();
223 }
224
225 build(dof_handler, external_partitioners);
226 }
227 else
228 {
229 build(dof_handler);
230 }
231}
232
233
234
235template <int dim, typename Number>
236void
238 const unsigned int to_level,
241{
242 dst = 0;
243 prolongate_and_add(to_level, dst, src);
244}
245
246
247
248template <int dim, typename Number>
249void
251 const unsigned int to_level,
254{
255 Assert((to_level >= 1) && (to_level <= level_dof_indices.size()),
256 ExcIndexRange(to_level, 1, level_dof_indices.size() + 1));
257
258 const bool src_inplace = src.get_partitioner().get() ==
259 this->vector_partitioners[to_level - 1].get();
260 if (src_inplace == false)
261 {
262 if (this->ghosted_level_vector[to_level - 1].get_partitioner().get() !=
263 this->vector_partitioners[to_level - 1].get())
264 this->ghosted_level_vector[to_level - 1].reinit(
265 this->vector_partitioners[to_level - 1]);
266 this->ghosted_level_vector[to_level - 1].copy_locally_owned_data_from(
267 src);
268 }
269
270 const bool dst_inplace =
271 dst.get_partitioner().get() == this->vector_partitioners[to_level].get();
272 if (dst_inplace == false)
273 {
274 if (this->ghosted_level_vector[to_level].get_partitioner().get() !=
275 this->vector_partitioners[to_level].get())
276 this->ghosted_level_vector[to_level].reinit(
277 this->vector_partitioners[to_level]);
278 AssertDimension(this->ghosted_level_vector[to_level].locally_owned_size(),
279 dst.locally_owned_size());
280 this->ghosted_level_vector[to_level] = 0.;
281 }
282
284 src_inplace ? src : this->ghosted_level_vector[to_level - 1];
286 dst_inplace ? dst : this->ghosted_level_vector[to_level];
287
288 src_vec.update_ghost_values();
289 // the implementation in do_prolongate_add is templated in the degree of the
290 // element (for efficiency reasons), so we need to find the appropriate
291 // kernel here...
292 if (fe_degree == 0)
293 do_prolongate_add<0>(to_level, dst_vec, src_vec);
294 else if (fe_degree == 1)
295 do_prolongate_add<1>(to_level, dst_vec, src_vec);
296 else if (fe_degree == 2)
297 do_prolongate_add<2>(to_level, dst_vec, src_vec);
298 else if (fe_degree == 3)
299 do_prolongate_add<3>(to_level, dst_vec, src_vec);
300 else if (fe_degree == 4)
301 do_prolongate_add<4>(to_level, dst_vec, src_vec);
302 else if (fe_degree == 5)
303 do_prolongate_add<5>(to_level, dst_vec, src_vec);
304 else if (fe_degree == 6)
305 do_prolongate_add<6>(to_level, dst_vec, src_vec);
306 else if (fe_degree == 7)
307 do_prolongate_add<7>(to_level, dst_vec, src_vec);
308 else if (fe_degree == 8)
309 do_prolongate_add<8>(to_level, dst_vec, src_vec);
310 else if (fe_degree == 9)
311 do_prolongate_add<9>(to_level, dst_vec, src_vec);
312 else if (fe_degree == 10)
313 do_prolongate_add<10>(to_level, dst_vec, src_vec);
314 else
315 do_prolongate_add<-1>(to_level, dst_vec, src_vec);
316
318 if (dst_inplace == false)
319 dst += dst_vec;
320
321 if (src_inplace == true)
323}
324
325
326
327template <int dim, typename Number>
328void
330 const unsigned int from_level,
333{
334 Assert((from_level >= 1) && (from_level <= level_dof_indices.size()),
335 ExcIndexRange(from_level, 1, level_dof_indices.size() + 1));
336
337 const bool src_inplace =
338 src.get_partitioner().get() == this->vector_partitioners[from_level].get();
339 if (src_inplace == false)
340 {
341 if (this->ghosted_level_vector[from_level].get_partitioner().get() !=
342 this->vector_partitioners[from_level].get())
343 this->ghosted_level_vector[from_level].reinit(
344 this->vector_partitioners[from_level]);
345 this->ghosted_level_vector[from_level].copy_locally_owned_data_from(src);
346 }
347
348 const bool dst_inplace = dst.get_partitioner().get() ==
349 this->vector_partitioners[from_level - 1].get();
350 if (dst_inplace == false)
351 {
352 if (this->ghosted_level_vector[from_level - 1].get_partitioner().get() !=
353 this->vector_partitioners[from_level - 1].get())
354 this->ghosted_level_vector[from_level - 1].reinit(
355 this->vector_partitioners[from_level - 1]);
357 this->ghosted_level_vector[from_level - 1].locally_owned_size(),
358 dst.locally_owned_size());
359 this->ghosted_level_vector[from_level - 1] = 0.;
360 }
361
363 src_inplace ? src : this->ghosted_level_vector[from_level];
365 dst_inplace ? dst : this->ghosted_level_vector[from_level - 1];
366
367 src_vec.update_ghost_values();
368
369 if (fe_degree == 0)
370 do_restrict_add<0>(from_level, dst_vec, src_vec);
371 else if (fe_degree == 1)
372 do_restrict_add<1>(from_level, dst_vec, src_vec);
373 else if (fe_degree == 2)
374 do_restrict_add<2>(from_level, dst_vec, src_vec);
375 else if (fe_degree == 3)
376 do_restrict_add<3>(from_level, dst_vec, src_vec);
377 else if (fe_degree == 4)
378 do_restrict_add<4>(from_level, dst_vec, src_vec);
379 else if (fe_degree == 5)
380 do_restrict_add<5>(from_level, dst_vec, src_vec);
381 else if (fe_degree == 6)
382 do_restrict_add<6>(from_level, dst_vec, src_vec);
383 else if (fe_degree == 7)
384 do_restrict_add<7>(from_level, dst_vec, src_vec);
385 else if (fe_degree == 8)
386 do_restrict_add<8>(from_level, dst_vec, src_vec);
387 else if (fe_degree == 9)
388 do_restrict_add<9>(from_level, dst_vec, src_vec);
389 else if (fe_degree == 10)
390 do_restrict_add<10>(from_level, dst_vec, src_vec);
391 else
392 // go to the non-templated version of the evaluator
393 do_restrict_add<-1>(from_level, dst_vec, src_vec);
394
396 if (dst_inplace == false)
397 dst += dst_vec;
398
399 if (src_inplace == true)
401}
402
403
404
405template <int dim, typename Number>
406template <int degree>
407void
409 const unsigned int to_level,
412{
413 const unsigned int vec_size = VectorizedArray<Number>::size();
414 const unsigned int degree_size = (degree > -1 ? degree : fe_degree) + 1;
415 const unsigned int n_child_dofs_1d = 2 * degree_size - element_is_continuous;
416 const unsigned int n_scalar_cell_dofs =
417 Utilities::fixed_power<dim>(n_child_dofs_1d);
418 constexpr unsigned int three_to_dim = Utilities::pow(3, dim);
419
420 // If we have user defined MG constraints, we must create
421 // a non-const, ghosted version of the source vector to distribute
422 // constraints.
425 if (this->mg_constrained_dofs != nullptr &&
426 this->mg_constrained_dofs->get_user_constraint_matrix(to_level - 1)
427 .get_local_lines()
428 .size() > 0)
429 {
431
432 // Distribute any user defined constraints
433 this->mg_constrained_dofs->get_user_constraint_matrix(to_level - 1)
434 .distribute(copy_src);
435
436 // Re-initialize new ghosted vector with correct constraints
437 new_src.reinit(copy_src);
438 new_src = copy_src;
439 new_src.update_ghost_values();
440
441 to_use = &new_src;
442 }
443
444 for (unsigned int cell = 0; cell < n_owned_level_cells[to_level - 1];
445 cell += vec_size)
446 {
447 const unsigned int n_lanes =
448 (cell + vec_size > n_owned_level_cells[to_level - 1]) ?
449 (n_owned_level_cells[to_level - 1] - cell) :
450 vec_size;
451
452 // read from source vector
453 for (unsigned int v = 0; v < n_lanes; ++v)
454 {
455 const unsigned int shift =
456 internal::MGTransfer::compute_shift_within_children<dim>(
457 parent_child_connect[to_level - 1][cell + v].second,
458 fe_degree + 1 - element_is_continuous,
459 fe_degree);
460 const unsigned int *indices =
461 &level_dof_indices[to_level - 1]
462 [parent_child_connect[to_level - 1][cell + v]
463 .first *
464 n_child_cell_dofs +
465 shift];
466 for (unsigned int c = 0, m = 0; c < n_components; ++c)
467 {
468 for (unsigned int k = 0; k < (dim > 2 ? degree_size : 1); ++k)
469 for (unsigned int j = 0; j < (dim > 1 ? degree_size : 1); ++j)
470 for (unsigned int i = 0; i < degree_size; ++i, ++m)
471 evaluation_data[m][v] = to_use->local_element(
472 indices[c * n_scalar_cell_dofs +
473 k * n_child_dofs_1d * n_child_dofs_1d +
474 j * n_child_dofs_1d + i]);
475
476 // apply Dirichlet boundary conditions on parent cell
477 for (std::vector<unsigned short>::const_iterator i =
478 dirichlet_indices[to_level - 1][cell + v].begin();
479 i != dirichlet_indices[to_level - 1][cell + v].end();
480 ++i)
481 evaluation_data[*i][v] = 0.;
482 }
483 }
484
485 AssertDimension(prolongation_matrix_1d.size(),
486 degree_size * n_child_dofs_1d);
487 // perform tensorized operation
488 if (element_is_continuous)
489 {
490 // must go through the components backwards because we want to write
491 // the output to the same array as the input
492 for (int c = n_components - 1; c >= 0; --c)
496 dim,
497 degree + 1,
498 2 * degree + 1>::do_forward(1,
499 prolongation_matrix_1d,
500 evaluation_data.begin() +
501 c * Utilities::fixed_power<dim>(
502 degree_size),
503 evaluation_data.begin() +
504 c * n_scalar_cell_dofs,
505 fe_degree + 1,
506 2 * fe_degree + 1);
508 degree != -1 ? 2 * degree + 1 :
509 -1,
511 &weights_on_refined[to_level - 1][(cell / vec_size) * three_to_dim],
512 n_components,
513 2 * fe_degree + 1,
514 evaluation_data.begin());
515 }
516 else
517 {
518 for (int c = n_components - 1; c >= 0; --c)
522 dim,
523 degree + 1,
524 2 * degree + 2>::do_forward(1,
525 prolongation_matrix_1d,
526 evaluation_data.begin() +
527 c * Utilities::fixed_power<dim>(
528 degree_size),
529 evaluation_data.begin() +
530 c * n_scalar_cell_dofs,
531 fe_degree + 1,
532 2 * fe_degree + 2);
533 }
534
535 // write into dst vector
536 const unsigned int *indices =
537 &level_dof_indices[to_level][cell * n_child_cell_dofs];
538 for (unsigned int v = 0; v < n_lanes; ++v)
539 {
540 for (unsigned int i = 0; i < n_child_cell_dofs; ++i)
541 dst.local_element(indices[i]) += evaluation_data[i][v];
542 indices += n_child_cell_dofs;
543 }
544 }
545}
546
547
548
549template <int dim, typename Number>
550template <int degree>
551void
553 const unsigned int from_level,
556{
557 const unsigned int vec_size = VectorizedArray<Number>::size();
558 const unsigned int degree_size = (degree > -1 ? degree : fe_degree) + 1;
559 const unsigned int n_child_dofs_1d = 2 * degree_size - element_is_continuous;
560 const unsigned int n_scalar_cell_dofs =
561 Utilities::fixed_power<dim>(n_child_dofs_1d);
562 constexpr unsigned int three_to_dim = Utilities::pow(3, dim);
563
564 for (unsigned int cell = 0; cell < n_owned_level_cells[from_level - 1];
565 cell += vec_size)
566 {
567 const unsigned int n_lanes =
568 (cell + vec_size > n_owned_level_cells[from_level - 1]) ?
569 (n_owned_level_cells[from_level - 1] - cell) :
570 vec_size;
571
572 // read from source vector
573 {
574 const unsigned int *indices =
575 &level_dof_indices[from_level][cell * n_child_cell_dofs];
576 for (unsigned int v = 0; v < n_lanes; ++v)
577 {
578 for (unsigned int i = 0; i < n_child_cell_dofs; ++i)
579 evaluation_data[i][v] = src.local_element(indices[i]);
580 indices += n_child_cell_dofs;
581 }
582 }
583
584 AssertDimension(prolongation_matrix_1d.size(),
585 degree_size * n_child_dofs_1d);
586 // perform tensorized operation
587 if (element_is_continuous)
588 {
590 degree != -1 ? 2 * degree + 1 :
591 -1,
593 &weights_on_refined[from_level - 1]
594 [(cell / vec_size) * three_to_dim],
595 n_components,
596 2 * fe_degree + 1,
597 evaluation_data.data());
598 for (unsigned int c = 0; c < n_components; ++c)
602 dim,
603 degree + 1,
604 2 * degree + 1>::do_backward(1,
605 prolongation_matrix_1d,
606 false,
607 evaluation_data.begin() +
608 c * n_scalar_cell_dofs,
609 evaluation_data.begin() +
610 c * Utilities::fixed_power<dim>(
611 degree_size),
612 fe_degree + 1,
613 2 * fe_degree + 1);
614 }
615 else
616 {
617 for (unsigned int c = 0; c < n_components; ++c)
621 dim,
622 degree + 1,
623 2 * degree + 2>::do_backward(1,
624 prolongation_matrix_1d,
625 false,
626 evaluation_data.begin() +
627 c * n_scalar_cell_dofs,
628 evaluation_data.begin() +
629 c * Utilities::fixed_power<dim>(
630 degree_size),
631 fe_degree + 1,
632 2 * fe_degree + 2);
633 }
634
635 // write into dst vector
636 for (unsigned int v = 0; v < n_lanes; ++v)
637 {
638 const unsigned int shift =
639 internal::MGTransfer::compute_shift_within_children<dim>(
640 parent_child_connect[from_level - 1][cell + v].second,
641 fe_degree + 1 - element_is_continuous,
642 fe_degree);
644 parent_child_connect[from_level - 1][cell + v].first *
645 n_child_cell_dofs +
646 n_child_cell_dofs - 1,
647 level_dof_indices[from_level - 1].size());
648 const unsigned int *indices =
649 &level_dof_indices[from_level - 1]
650 [parent_child_connect[from_level - 1][cell + v]
651 .first *
652 n_child_cell_dofs +
653 shift];
654 for (unsigned int c = 0, m = 0; c < n_components; ++c)
655 {
656 // apply Dirichlet boundary conditions on parent cell
657 for (std::vector<unsigned short>::const_iterator i =
658 dirichlet_indices[from_level - 1][cell + v].begin();
659 i != dirichlet_indices[from_level - 1][cell + v].end();
660 ++i)
661 evaluation_data[*i][v] = 0.;
662
663 for (unsigned int k = 0; k < (dim > 2 ? degree_size : 1); ++k)
664 for (unsigned int j = 0; j < (dim > 1 ? degree_size : 1); ++j)
665 for (unsigned int i = 0; i < degree_size; ++i, ++m)
666 dst.local_element(
667 indices[c * n_scalar_cell_dofs +
668 k * n_child_dofs_1d * n_child_dofs_1d +
669 j * n_child_dofs_1d + i]) +=
670 evaluation_data[m][v];
671 }
672 }
673 }
674}
675
676
677
678template <int dim, typename Number>
679std::size_t
681{
682 std::size_t memory = MGLevelGlobalTransfer<
683 LinearAlgebra::distributed::Vector<Number>>::memory_consumption();
684 memory += MemoryConsumption::memory_consumption(level_dof_indices);
685 memory += MemoryConsumption::memory_consumption(parent_child_connect);
686 memory += MemoryConsumption::memory_consumption(n_owned_level_cells);
687 memory += MemoryConsumption::memory_consumption(prolongation_matrix_1d);
688 memory += MemoryConsumption::memory_consumption(evaluation_data);
689 memory += MemoryConsumption::memory_consumption(weights_on_refined);
690 memory += MemoryConsumption::memory_consumption(dirichlet_indices);
691 return memory;
692}
693
694
695
696template <int dim, typename Number>
698 const MGConstrainedDoFs &mg_c)
700 Number,
701 MGTransferMatrixFree<dim, Number>>(true)
702{
703 this->matrix_free_transfer_vector.emplace_back(mg_c);
704}
705
706
707
708template <int dim, typename Number>
710 const std::vector<MGConstrainedDoFs> &mg_c)
712 Number,
713 MGTransferMatrixFree<dim, Number>>(false)
714{
715 for (const auto &constrained_block_dofs : mg_c)
716 this->matrix_free_transfer_vector.emplace_back(constrained_block_dofs);
717}
718
719
720
721template <int dim, typename Number>
722void
724 const MGConstrainedDoFs &mg_c)
725{
726 Assert(this->same_for_all,
727 ExcMessage("This object was initialized with support for usage with "
728 "one DoFHandler for each block, but this method assumes "
729 "that the same DoFHandler is used for all the blocks!"));
730 AssertDimension(this->matrix_free_transfer_vector.size(), 1);
731
732 this->matrix_free_transfer_vector[0].initialize_constraints(mg_c);
733}
734
735
736
737template <int dim, typename Number>
738void
740 const std::vector<MGConstrainedDoFs> &mg_c)
741{
742 Assert(!this->same_for_all,
743 ExcMessage("This object was initialized with support for using "
744 "the same DoFHandler for all the blocks, but this "
745 "method assumes that there is a separate DoFHandler "
746 "for each block!"));
747 AssertDimension(this->matrix_free_transfer_vector.size(), mg_c.size());
748
749 for (unsigned int i = 0; i < mg_c.size(); ++i)
750 this->matrix_free_transfer_vector[i].initialize_constraints(mg_c[i]);
751}
752
753
754
755template <int dim, typename Number>
756void
758 const DoFHandler<dim> &dof_handler)
759{
760 AssertDimension(this->matrix_free_transfer_vector.size(), 1);
761 this->matrix_free_transfer_vector[0].build(dof_handler);
762}
763
764
765
766template <int dim, typename Number>
767void
769 const std::vector<const DoFHandler<dim> *> &dof_handler)
770{
771 AssertDimension(this->matrix_free_transfer_vector.size(), dof_handler.size());
772 for (unsigned int i = 0; i < dof_handler.size(); ++i)
773 this->matrix_free_transfer_vector[i].build(*dof_handler[i]);
774}
775
776
777
778template <int dim, typename Number>
779std::size_t
781{
782 std::size_t total_memory_consumption = 0;
783 for (const auto &el : this->matrix_free_transfer_vector)
784 total_memory_consumption += el.memory_consumption();
785 return total_memory_consumption;
786}
787
788
789
790template <int dim, typename Number>
793 const unsigned int b) const
794{
795 return matrix_free_transfer_vector[b];
796}
797
798
799
800template <int dim, typename Number>
801void
803{
804 this->matrix_free_transfer_vector.clear();
805}
806
807
808
809// explicit instantiation
810#include "mg_transfer_matrix_free.inst"
811
812
const Triangulation< dim, spacedim > & get_triangulation() const
bool has_level_dofs() const
Number local_element(const size_type local_index) const
size_type locally_owned_size() const
void compress(VectorOperation::values operation)
const std::shared_ptr< const Utilities::MPI::Partitioner > & get_partitioner() const
void reinit(const size_type size, const bool omit_zeroing_entries=false)
ObserverPointer< const MGConstrainedDoFs > mg_constrained_dofs
std::vector< MGTransferMatrixFree< dim, Number > > matrix_free_transfer_vector
void build(const DoFHandler< dim > &dof_handler)
MGTransferBlockMatrixFree()=default
const MGTransferMatrixFree< dim, Number > & get_matrix_free_transfer(const unsigned int b) const override
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
void do_restrict_add(const unsigned int from_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
void do_prolongate_add(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
virtual void prolongate(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const override
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
void build(const DoFHandler< dim > &dof_handler, const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > &external_partitioners=std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > >())
virtual void prolongate_and_add(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const override
std::size_t memory_consumption() const
virtual void restrict_and_add(const unsigned int from_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const override
virtual unsigned int n_global_levels() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
Point< 2 > second
Definition grid_out.cc:4624
Point< 2 > first
Definition grid_out.cc:4623
unsigned int level
Definition grid_out.cc:4626
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcIndexRange(std::size_t arg1, std::size_t arg2, std::size_t arg3)
static ::ExceptionBase & ExcMessage(std::string arg1)
std::size_t size
Definition mpi.cc:734
types::global_dof_index locally_owned_size
Definition mpi.cc:822
std::enable_if_t< std::is_fundamental_v< T >, std::size_t > memory_consumption(const T &t)
constexpr T pow(const T base, const int iexp)
Definition utilities.h:966
void weight_fe_q_dofs_by_entity(const Number *weights, const unsigned int n_components, const int n_points_1d_non_template, Number *data)