From: Martin Kronbichler Date: Thu, 19 Nov 2020 19:05:49 +0000 (+0100) Subject: Fix tests with MGTransferMatrixFree::interpolate_to_mg X-Git-Tag: v9.3.0-rc1~882^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=48f5be472267028d7be511380e911b9ffae49f06;p=dealii.git Fix tests with MGTransferMatrixFree::interpolate_to_mg --- diff --git a/include/deal.II/multigrid/mg_transfer_matrix_free.h b/include/deal.II/multigrid/mg_transfer_matrix_free.h index 7e86e2af75..48cdb2e5e4 100644 --- a/include/deal.II/multigrid/mg_transfer_matrix_free.h +++ b/include/deal.II/multigrid/mg_transfer_matrix_free.h @@ -505,6 +505,12 @@ MGTransferMatrixFree::interpolate_to_mg( const FiniteElement &fe = dof_handler.get_fe(); + for (unsigned int level = min_level; level <= max_level; ++level) + if (dst[level].size() != dof_handler.n_dofs(level) || + dst[level].local_size() != + dof_handler.locally_owned_mg_dofs(level).n_elements()) + dst[level].reinit(this->vector_partitioners[level]); + // copy fine level vector to active cells in MG hierarchy this->copy_to_mg(dof_handler, dst, src, true); @@ -513,15 +519,22 @@ MGTransferMatrixFree::interpolate_to_mg( // constraints... // do the transfer from level to level-1: + dst[max_level].update_ghost_values(); for (unsigned int level = max_level; level > min_level; --level) { // auxiliary vector which always has ghost elements - LinearAlgebra::distributed::Vector ghosted_fine( - this->vector_partitioners[level]); - ghosted_fine.copy_locally_owned_data_from(dst[level]); - ghosted_fine.update_ghost_values(); - LinearAlgebra::distributed::Vector ghosted_coarse( - this->vector_partitioners[level - 1]); + const LinearAlgebra::distributed::Vector *input = nullptr; + LinearAlgebra::distributed::Vector ghosted_fine; + if (dst[level].get_partitioner().get() == + this->vector_partitioners[level].get()) + input = &dst[level]; + else + { + ghosted_fine.reinit(this->vector_partitioners[level]); + ghosted_fine.copy_locally_owned_data_from(dst[level]); + ghosted_fine.update_ghost_values(); + input = &ghosted_fine; + } std::vector dof_values_coarse(fe.n_dofs_per_cell()); Vector dof_values_fine(fe.n_dofs_per_cell()); @@ -541,7 +554,7 @@ MGTransferMatrixFree::interpolate_to_mg( { cell->child(child)->get_mg_dof_indices(dof_indices); for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i) - dof_values_fine(i) = ghosted_fine(dof_indices[i]); + dof_values_fine(i) = (*input)(dof_indices[i]); fe.get_restriction_matrix(child, cell->refinement_case()) .vmult(tmp, dof_values_fine); for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i) @@ -552,9 +565,12 @@ MGTransferMatrixFree::interpolate_to_mg( } cell->get_mg_dof_indices(dof_indices); for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i) - ghosted_coarse(dof_indices[i]) = dof_values_coarse[i]; + if (dof_handler.locally_owned_mg_dofs(level - 1).is_element( + dof_indices[i])) + dst[level - 1](dof_indices[i]) = dof_values_coarse[i]; } - dst[level - 1].copy_locally_owned_data_from(ghosted_coarse); + + dst[level - 1].update_ghost_values(); } } diff --git a/tests/mappings/mapping_q_eulerian_08.cc b/tests/mappings/mapping_q_eulerian_08.cc index 0fc87e44fb..209c064f5b 100644 --- a/tests/mappings/mapping_q_eulerian_08.cc +++ b/tests/mappings/mapping_q_eulerian_08.cc @@ -199,7 +199,7 @@ test(const unsigned int n_ref = 0) MGTransferMatrixFree mg_transfer_euler( mg_constrained_dofs_euler); - mg_transfer_euler.build(dof_handler); + mg_transfer_euler.build(dof_handler_euler); // now the core of the test: const unsigned int max_level = diff --git a/tests/matrix_free/interpolate_to_mg.cc b/tests/matrix_free/interpolate_to_mg.cc index 07cdc064ba..f09b570b0c 100644 --- a/tests/matrix_free/interpolate_to_mg.cc +++ b/tests/matrix_free/interpolate_to_mg.cc @@ -151,7 +151,7 @@ test(const unsigned int n_glob_ref = 2, const unsigned int n_ref = 0) fine_projection.update_ghost_values(); // output for debug purposes: - if (true) + if (false) { DataOut data_out; data_out.attach_dof_handler(dof_handler); @@ -190,6 +190,14 @@ test(const unsigned int n_glob_ref = 2, const unsigned int n_ref = 0) const unsigned int min_level = 0; MGLevelObject> level_projection(min_level, max_level); + for (unsigned int level = min_level; level <= max_level; ++level) + { + IndexSet set; + DoFTools::extract_locally_relevant_level_dofs(dof_handler, level, set); + level_projection[level].reinit(dof_handler.locally_owned_mg_dofs(level), + set, + mpi_communicator); + } mg_transfer.interpolate_to_mg(dof_handler, level_projection, fine_projection); // now go through all GMG levels and make sure FE field can represent