36#include <deal.II/multigrid/mg_transfer_matrix_free.templates.h>
43template <
int dim,
typename Number>
46 , element_is_continuous(false)
48 , n_child_cell_dofs(0)
53template <
int dim,
typename Number>
57 , element_is_continuous(false)
59 , n_child_cell_dofs(0)
66template <
int dim,
typename Number>
71 this->mg_constrained_dofs = &mg_c;
76template <
int dim,
typename Number>
83 element_is_continuous =
false;
85 n_child_cell_dofs = 0;
86 level_dof_indices.clear();
87 parent_child_connect.clear();
88 dirichlet_indices.clear();
89 n_owned_level_cells.clear();
90 prolongation_matrix_1d.clear();
91 evaluation_data.clear();
92 weights_on_refined.clear();
97template <
int dim,
typename Number>
101 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
102 &external_partitioners)
106 "The underlying DoFHandler object has not had its "
107 "distribute_mg_dofs() function called, but this is a prerequisite "
108 "for multigrid transfers. You will need to call this function, "
109 "probably close to where you already call distribute_dofs()."));
110 if (external_partitioners.size() > 0)
112 Assert(this->initialize_dof_vector ==
nullptr,
113 ExcMessage(
"An initialize_dof_vector function has already "
114 "been registered in the constructor!"));
116 this->initialize_dof_vector =
117 [external_partitioners](
118 const unsigned int level,
120 vec.reinit(external_partitioners[
level]);
124 this->fill_and_communicate_copy_indices(dof_handler);
126 vector_partitioners.resize(0,
129 for (
unsigned int level = 0;
level <= this->ghosted_level_vector.max_level();
131 vector_partitioners[
level] =
132 this->ghosted_level_vector[
level].get_partitioner();
134 std::vector<std::vector<Number>> weights_unvectorized;
138 internal::MGTransfer::setup_transfer<dim, Number>(
140 this->mg_constrained_dofs,
141 external_partitioners,
144 parent_child_connect,
147 weights_unvectorized,
148 this->copy_indices_global_mine,
149 vector_partitioners);
153 this->ghosted_level_vector.resize(0, vector_partitioners.max_level());
154 for (
unsigned int level = 0;
level <= vector_partitioners.max_level();
156 if (external_partitioners.size() == vector_partitioners.max_level() + 1 &&
157 external_partitioners[
level].get() == vector_partitioners[
level].get())
158 this->ghosted_level_vector[
level].reinit(0);
160 this->ghosted_level_vector[
level].reinit(vector_partitioners[
level]);
175 const unsigned int n_levels =
178 const unsigned int n_weights_per_cell = Utilities::fixed_power<dim>(3);
179 weights_on_refined.resize(n_levels - 1);
182 weights_on_refined[
level - 1].resize(
183 ((n_owned_level_cells[
level - 1] + vec_size - 1) / vec_size) *
186 for (
unsigned int c = 0; c < n_owned_level_cells[
level - 1]; ++c)
188 const unsigned int comp = c / vec_size;
189 const unsigned int v = c % vec_size;
190 for (
unsigned int i = 0; i < n_weights_per_cell; ++i)
192 weights_on_refined[
level - 1][comp * n_weights_per_cell + i][v] =
193 weights_unvectorized[
level - 1][c * n_weights_per_cell + i];
198 evaluation_data.resize(n_child_cell_dofs);
203template <
int dim,
typename Number>
207 const std::function<
void(
const unsigned int,
209 &initialize_dof_vector)
211 if (initialize_dof_vector)
213 const unsigned int n_levels =
216 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
217 external_partitioners(n_levels);
222 initialize_dof_vector(
level, vector);
226 build(dof_handler, external_partitioners);
236template <
int dim,
typename Number>
239 const unsigned int to_level,
244 prolongate_and_add(to_level, dst, src);
249template <
int dim,
typename Number>
252 const unsigned int to_level,
256 Assert((to_level >= 1) && (to_level <= level_dof_indices.size()),
260 this->vector_partitioners[to_level - 1].get();
261 if (src_inplace ==
false)
263 if (this->ghosted_level_vector[to_level - 1].get_partitioner().get() !=
264 this->vector_partitioners[to_level - 1].get())
265 this->ghosted_level_vector[to_level - 1].reinit(
266 this->vector_partitioners[to_level - 1]);
267 this->ghosted_level_vector[to_level - 1].copy_locally_owned_data_from(
271 const bool dst_inplace =
272 dst.
get_partitioner().get() == this->vector_partitioners[to_level].get();
273 if (dst_inplace ==
false)
275 if (this->ghosted_level_vector[to_level].get_partitioner().get() !=
276 this->vector_partitioners[to_level].get())
277 this->ghosted_level_vector[to_level].reinit(
278 this->vector_partitioners[to_level]);
279 AssertDimension(this->ghosted_level_vector[to_level].locally_owned_size(),
281 this->ghosted_level_vector[to_level] = 0.;
285 src_inplace ? src : this->ghosted_level_vector[to_level - 1];
287 dst_inplace ? dst : this->ghosted_level_vector[to_level];
294 do_prolongate_add<0>(to_level, dst_vec, src_vec);
295 else if (fe_degree == 1)
296 do_prolongate_add<1>(to_level, dst_vec, src_vec);
297 else if (fe_degree == 2)
298 do_prolongate_add<2>(to_level, dst_vec, src_vec);
299 else if (fe_degree == 3)
300 do_prolongate_add<3>(to_level, dst_vec, src_vec);
301 else if (fe_degree == 4)
302 do_prolongate_add<4>(to_level, dst_vec, src_vec);
303 else if (fe_degree == 5)
304 do_prolongate_add<5>(to_level, dst_vec, src_vec);
305 else if (fe_degree == 6)
306 do_prolongate_add<6>(to_level, dst_vec, src_vec);
307 else if (fe_degree == 7)
308 do_prolongate_add<7>(to_level, dst_vec, src_vec);
309 else if (fe_degree == 8)
310 do_prolongate_add<8>(to_level, dst_vec, src_vec);
311 else if (fe_degree == 9)
312 do_prolongate_add<9>(to_level, dst_vec, src_vec);
313 else if (fe_degree == 10)
314 do_prolongate_add<10>(to_level, dst_vec, src_vec);
316 do_prolongate_add<-1>(to_level, dst_vec, src_vec);
319 if (dst_inplace ==
false)
322 if (src_inplace ==
true)
328template <
int dim,
typename Number>
331 const unsigned int from_level,
335 Assert((from_level >= 1) && (from_level <= level_dof_indices.size()),
338 const bool src_inplace =
339 src.
get_partitioner().get() == this->vector_partitioners[from_level].get();
340 if (src_inplace ==
false)
342 if (this->ghosted_level_vector[from_level].get_partitioner().get() !=
343 this->vector_partitioners[from_level].get())
344 this->ghosted_level_vector[from_level].reinit(
345 this->vector_partitioners[from_level]);
346 this->ghosted_level_vector[from_level].copy_locally_owned_data_from(src);
350 this->vector_partitioners[from_level - 1].get();
351 if (dst_inplace ==
false)
353 if (this->ghosted_level_vector[from_level - 1].get_partitioner().get() !=
354 this->vector_partitioners[from_level - 1].get())
355 this->ghosted_level_vector[from_level - 1].reinit(
356 this->vector_partitioners[from_level - 1]);
358 this->ghosted_level_vector[from_level - 1].locally_owned_size(),
360 this->ghosted_level_vector[from_level - 1] = 0.;
364 src_inplace ? src : this->ghosted_level_vector[from_level];
366 dst_inplace ? dst : this->ghosted_level_vector[from_level - 1];
371 do_restrict_add<0>(from_level, dst_vec, src_vec);
372 else if (fe_degree == 1)
373 do_restrict_add<1>(from_level, dst_vec, src_vec);
374 else if (fe_degree == 2)
375 do_restrict_add<2>(from_level, dst_vec, src_vec);
376 else if (fe_degree == 3)
377 do_restrict_add<3>(from_level, dst_vec, src_vec);
378 else if (fe_degree == 4)
379 do_restrict_add<4>(from_level, dst_vec, src_vec);
380 else if (fe_degree == 5)
381 do_restrict_add<5>(from_level, dst_vec, src_vec);
382 else if (fe_degree == 6)
383 do_restrict_add<6>(from_level, dst_vec, src_vec);
384 else if (fe_degree == 7)
385 do_restrict_add<7>(from_level, dst_vec, src_vec);
386 else if (fe_degree == 8)
387 do_restrict_add<8>(from_level, dst_vec, src_vec);
388 else if (fe_degree == 9)
389 do_restrict_add<9>(from_level, dst_vec, src_vec);
390 else if (fe_degree == 10)
391 do_restrict_add<10>(from_level, dst_vec, src_vec);
394 do_restrict_add<-1>(from_level, dst_vec, src_vec);
397 if (dst_inplace ==
false)
400 if (src_inplace ==
true)
406template <
int dim,
typename Number>
410 const unsigned int to_level,
415 const unsigned int degree_size = (degree > -1 ? degree : fe_degree) + 1;
416 const unsigned int n_child_dofs_1d = 2 * degree_size - element_is_continuous;
417 const unsigned int n_scalar_cell_dofs =
418 Utilities::fixed_power<dim>(n_child_dofs_1d);
426 if (this->mg_constrained_dofs !=
nullptr &&
427 this->mg_constrained_dofs->get_user_constraint_matrix(to_level - 1)
434 this->mg_constrained_dofs->get_user_constraint_matrix(to_level - 1)
435 .distribute(copy_src);
445 for (
unsigned int cell = 0; cell < n_owned_level_cells[to_level - 1];
448 const unsigned int n_lanes =
449 (cell + vec_size > n_owned_level_cells[to_level - 1]) ?
450 (n_owned_level_cells[to_level - 1] - cell) :
454 for (
unsigned int v = 0; v < n_lanes; ++v)
456 const unsigned int shift =
457 internal::MGTransfer::compute_shift_within_children<dim>(
458 parent_child_connect[to_level - 1][cell + v].
second,
459 fe_degree + 1 - element_is_continuous,
461 const unsigned int *indices =
462 &level_dof_indices[to_level - 1]
463 [parent_child_connect[to_level - 1][cell + v]
467 for (
unsigned int c = 0, m = 0; c < n_components; ++c)
469 for (
unsigned int k = 0; k < (dim > 2 ? degree_size : 1); ++k)
470 for (
unsigned int j = 0; j < (dim > 1 ? degree_size : 1); ++j)
471 for (
unsigned int i = 0; i < degree_size; ++i, ++m)
473 indices[c * n_scalar_cell_dofs +
474 k * n_child_dofs_1d * n_child_dofs_1d +
475 j * n_child_dofs_1d + i]);
478 for (std::vector<unsigned short>::const_iterator i =
479 dirichlet_indices[to_level - 1][cell + v].begin();
480 i != dirichlet_indices[to_level - 1][cell + v].end();
482 evaluation_data[*i][v] = 0.;
487 degree_size * n_child_dofs_1d);
489 if (element_is_continuous)
493 for (
int c = n_components - 1; c >= 0; --c)
499 2 * degree + 1>::do_forward(1,
500 prolongation_matrix_1d,
501 evaluation_data.begin() +
502 c * Utilities::fixed_power<dim>(
504 evaluation_data.begin() +
505 c * n_scalar_cell_dofs,
509 degree != -1 ? 2 * degree + 1 :
512 &weights_on_refined[to_level - 1][(cell / vec_size) * three_to_dim],
515 evaluation_data.begin());
519 for (
int c = n_components - 1; c >= 0; --c)
525 2 * degree + 2>::do_forward(1,
526 prolongation_matrix_1d,
527 evaluation_data.begin() +
528 c * Utilities::fixed_power<dim>(
530 evaluation_data.begin() +
531 c * n_scalar_cell_dofs,
537 const unsigned int *indices =
538 &level_dof_indices[to_level][cell * n_child_cell_dofs];
539 for (
unsigned int v = 0; v < n_lanes; ++v)
541 for (
unsigned int i = 0; i < n_child_cell_dofs; ++i)
543 indices += n_child_cell_dofs;
550template <
int dim,
typename Number>
554 const unsigned int from_level,
559 const unsigned int degree_size = (degree > -1 ? degree : fe_degree) + 1;
560 const unsigned int n_child_dofs_1d = 2 * degree_size - element_is_continuous;
561 const unsigned int n_scalar_cell_dofs =
562 Utilities::fixed_power<dim>(n_child_dofs_1d);
565 for (
unsigned int cell = 0; cell < n_owned_level_cells[from_level - 1];
568 const unsigned int n_lanes =
569 (cell + vec_size > n_owned_level_cells[from_level - 1]) ?
570 (n_owned_level_cells[from_level - 1] - cell) :
575 const unsigned int *indices =
576 &level_dof_indices[from_level][cell * n_child_cell_dofs];
577 for (
unsigned int v = 0; v < n_lanes; ++v)
579 for (
unsigned int i = 0; i < n_child_cell_dofs; ++i)
581 indices += n_child_cell_dofs;
586 degree_size * n_child_dofs_1d);
588 if (element_is_continuous)
591 degree != -1 ? 2 * degree + 1 :
594 &weights_on_refined[from_level - 1]
595 [(cell / vec_size) * three_to_dim],
598 evaluation_data.data());
599 for (
unsigned int c = 0; c < n_components; ++c)
605 2 * degree + 1>::do_backward(1,
606 prolongation_matrix_1d,
608 evaluation_data.begin() +
609 c * n_scalar_cell_dofs,
610 evaluation_data.begin() +
611 c * Utilities::fixed_power<dim>(
618 for (
unsigned int c = 0; c < n_components; ++c)
624 2 * degree + 2>::do_backward(1,
625 prolongation_matrix_1d,
627 evaluation_data.begin() +
628 c * n_scalar_cell_dofs,
629 evaluation_data.begin() +
630 c * Utilities::fixed_power<dim>(
637 for (
unsigned int v = 0; v < n_lanes; ++v)
639 const unsigned int shift =
640 internal::MGTransfer::compute_shift_within_children<dim>(
641 parent_child_connect[from_level - 1][cell + v].
second,
642 fe_degree + 1 - element_is_continuous,
645 parent_child_connect[from_level - 1][cell + v].
first *
647 n_child_cell_dofs - 1,
648 level_dof_indices[from_level - 1].size());
649 const unsigned int *indices =
650 &level_dof_indices[from_level - 1]
651 [parent_child_connect[from_level - 1][cell + v]
655 for (
unsigned int c = 0, m = 0; c < n_components; ++c)
658 for (std::vector<unsigned short>::const_iterator i =
659 dirichlet_indices[from_level - 1][cell + v].begin();
660 i != dirichlet_indices[from_level - 1][cell + v].end();
662 evaluation_data[*i][v] = 0.;
664 for (
unsigned int k = 0; k < (dim > 2 ? degree_size : 1); ++k)
665 for (
unsigned int j = 0; j < (dim > 1 ? degree_size : 1); ++j)
666 for (
unsigned int i = 0; i < degree_size; ++i, ++m)
668 indices[c * n_scalar_cell_dofs +
669 k * n_child_dofs_1d * n_child_dofs_1d +
670 j * n_child_dofs_1d + i]) +=
671 evaluation_data[m][v];
679template <
int dim,
typename Number>
697template <
int dim,
typename Number>
709template <
int dim,
typename Number>
711 const std::vector<MGConstrainedDoFs> &mg_c)
716 for (
const auto &constrained_block_dofs : mg_c)
722template <
int dim,
typename Number>
727 Assert(this->same_for_all,
728 ExcMessage(
"This object was initialized with support for usage with "
729 "one DoFHandler for each block, but this method assumes "
730 "that the same DoFHandler is used for all the blocks!"));
733 this->matrix_free_transfer_vector[0].initialize_constraints(mg_c);
738template <
int dim,
typename Number>
741 const std::vector<MGConstrainedDoFs> &mg_c)
743 Assert(!this->same_for_all,
744 ExcMessage(
"This object was initialized with support for using "
745 "the same DoFHandler for all the blocks, but this "
746 "method assumes that there is a separate DoFHandler "
748 AssertDimension(this->matrix_free_transfer_vector.size(), mg_c.size());
750 for (
unsigned int i = 0; i < mg_c.size(); ++i)
751 this->matrix_free_transfer_vector[i].initialize_constraints(mg_c[i]);
756template <
int dim,
typename Number>
762 this->matrix_free_transfer_vector[0].build(dof_handler);
767template <
int dim,
typename Number>
772 AssertDimension(this->matrix_free_transfer_vector.size(), dof_handler.size());
773 for (
unsigned int i = 0; i < dof_handler.size(); ++i)
774 this->matrix_free_transfer_vector[i].build(*dof_handler[i]);
779template <
int dim,
typename Number>
783 std::size_t total_memory_consumption = 0;
784 for (
const auto &el : this->matrix_free_transfer_vector)
785 total_memory_consumption += el.memory_consumption();
786 return total_memory_consumption;
791template <
int dim,
typename Number>
794 const unsigned int b)
const
796 return matrix_free_transfer_vector[b];
801template <
int dim,
typename Number>
805 this->matrix_free_transfer_vector.clear();
811#include "mg_transfer_matrix_free.inst"
const Triangulation< dim, spacedim > & get_triangulation() const
bool has_level_dofs() const
void zero_out_ghost_values() const
Number local_element(const size_type local_index) const
void update_ghost_values() const
size_type locally_owned_size() const
const std::shared_ptr< const Utilities::MPI::Partitioner > & get_partitioner() const
virtual void compress(VectorOperation::values operation) override
void reinit(const size_type size, const bool omit_zeroing_entries=false)
SmartPointer< const MGConstrainedDoFs, MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > > > mg_constrained_dofs
std::vector< MGTransferMatrixFree< dim, Number > > matrix_free_transfer_vector
std::size_t memory_consumption() const
MGTransferBlockMatrixFree()=default
const MGTransferMatrixFree< dim, Number > & get_matrix_free_transfer(const unsigned int b) const override
void build(const DoFHandler< dim, dim > &dof_handler)
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, 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
static constexpr std::size_t size()
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#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::enable_if_t< std::is_fundamental< T >::value, std::size_t > memory_consumption(const T &t)
constexpr T pow(const T base, const int iexp)
void weight_fe_q_dofs_by_entity(const Number *weights, const unsigned int n_components, const int n_points_1d_non_template, Number *data)
std::vector< Number > prolongation_matrix_1d
unsigned int n_child_cell_dofs
bool element_is_continuous
unsigned int n_components