36#include <deal.II/multigrid/mg_transfer.templates.h>
47template <
typename VectorType>
48template <
int dim,
int spacedim>
56 copy_indices_global_mine,
57 copy_indices_level_mine);
61 bool my_perform_plain_copy =
67 if (my_perform_plain_copy)
75 for (
const auto ©_index : copy_indices.back())
76 if (copy_index.first != copy_index.second)
78 my_perform_plain_copy =
false;
89 ptria->get_communicator()) == 1);
91 perform_plain_copy = my_perform_plain_copy;
96template <
typename VectorType>
101 copy_indices.clear();
102 copy_indices_global_mine.clear();
103 copy_indices_level_mine.clear();
104 component_to_block_map.resize(0);
105 mg_constrained_dofs =
nullptr;
106 perform_plain_copy =
false;
111template <
typename VectorType>
117 for (
unsigned int i = 0; i < copy_indices[
level].size(); ++i)
118 os <<
"copy_indices[" <<
level <<
"]\t" << copy_indices[
level][i].
first
119 <<
'\t' << copy_indices[
level][i].
second << std::endl;
122 for (
unsigned int level = 0;
level < copy_indices_level_mine.size(); ++
level)
124 for (
unsigned int i = 0; i < copy_indices_level_mine[
level].size(); ++i)
125 os <<
"copy_ifrom [" <<
level <<
"]\t"
126 << copy_indices_level_mine[
level][i].
first <<
'\t'
127 << copy_indices_level_mine[
level][i].
second << std::endl;
129 for (
unsigned int level = 0;
level < copy_indices_global_mine.size(); ++
level)
131 for (
unsigned int i = 0; i < copy_indices_global_mine[
level].size(); ++i)
132 os <<
"copy_ito [" <<
level <<
"]\t"
133 << copy_indices_global_mine[
level][i].
first <<
'\t'
134 << copy_indices_global_mine[
level][i].
second << std::endl;
140template <
typename VectorType>
144 std::size_t result =
sizeof(*this);
159 template <
int dim,
int spacedim,
typename Number>
165 const bool transfer_solution_vectors,
171 &ghosted_level_vector)
175 std::vector<std::pair<types::global_dof_index, types::global_dof_index>>>
178 std::vector<std::pair<types::global_dof_index, types::global_dof_index>>>
179 my_copy_indices_global_mine;
181 std::vector<std::pair<types::global_dof_index, types::global_dof_index>>>
182 my_copy_indices_level_mine;
187 my_copy_indices_global_mine,
188 my_copy_indices_level_mine,
189 !transfer_solution_vectors);
201 std::vector<types::global_dof_index> accessed_indices;
202 ghosted_level_vector.resize(0,
205 std::vector<IndexSet> level_index_set(
210 for (
const auto &indices : my_copy_indices_level_mine[l])
211 accessed_indices.push_back(indices.
first);
212 std::vector<types::global_dof_index> accessed_level_indices;
213 for (
const auto &indices : my_copy_indices_global_mine[l])
214 accessed_level_indices.push_back(indices.
second);
215 std::sort(accessed_level_indices.begin(), accessed_level_indices.end());
217 level_index_set[l].add_indices(accessed_level_indices.begin(),
218 accessed_level_indices.end());
219 level_index_set[l].compress();
224 std::sort(accessed_indices.begin(), accessed_indices.end());
225 index_set.add_indices(accessed_indices.begin(), accessed_indices.end());
226 index_set.compress();
235 copy_indices_level_mine.resize(
237 copy_indices_global_mine.resize(
239 for (
unsigned int level = 0;
246 *ghosted_level_vector[
level].get_partitioner();
248 auto translate_indices =
249 [&](
const std::vector<
250 std::pair<types::global_dof_index, types::global_dof_index>>
251 &global_copy_indices,
253 local_copy_indices.
reinit(2, global_copy_indices.size());
254 for (
unsigned int i = 0; i < global_copy_indices.size(); ++i)
257 global_copy_indices[i].
first);
259 global_copy_indices[i].
second);
264 translate_indices(my_copy_indices[
level], copy_indices[
level]);
267 translate_indices(my_copy_indices_level_mine[
level],
268 copy_indices_level_mine[
level]);
271 translate_indices(my_copy_indices_global_mine[
level],
272 copy_indices_global_mine[
level]);
277template <
typename Number>
278template <
int dim,
int spacedim>
285 fill_internal(mg_dof,
290 this->copy_indices_global_mine,
291 this->copy_indices_level_mine,
292 ghosted_global_vector,
293 ghosted_level_vector);
299 int have_refinement_edge_dofs = 0;
300 if (mg_constrained_dofs !=
nullptr)
301 for (
unsigned int level = 0;
304 if (mg_constrained_dofs->get_refinement_edge_indices(
level).n_elements() >
307 have_refinement_edge_dofs = 1;
313 std::vector<Table<2, unsigned int>> solution_copy_indices_global_mine;
315 solution_ghosted_level_vector;
317 fill_internal(mg_dof,
321 this->solution_copy_indices,
322 solution_copy_indices_global_mine,
323 this->solution_copy_indices_level_mine,
324 solution_ghosted_global_vector,
325 solution_ghosted_level_vector);
329 this->solution_copy_indices = this->copy_indices;
330 this->solution_copy_indices_level_mine = this->copy_indices_level_mine;
331 solution_ghosted_global_vector = ghosted_global_vector;
340 const bool my_perform_renumbered_plain_copy =
341 (this->copy_indices.back().n_cols() ==
343 (this->copy_indices_global_mine.back().n_rows() == 0) &&
344 (this->copy_indices_level_mine.back().n_rows() == 0);
346 bool my_perform_plain_copy =
false;
347 if (my_perform_renumbered_plain_copy)
349 my_perform_plain_copy =
true;
353 for (
unsigned int i = 0; i < this->copy_indices.back().n_cols(); ++i)
354 if (this->copy_indices.back()(0, i) != this->copy_indices.back()(1, i))
356 my_perform_plain_copy =
false;
366 perform_renumbered_plain_copy =
371 if (perform_renumbered_plain_copy)
373 for (
unsigned int i = 0; i < this->copy_indices.back().n_cols(); ++i)
376 ghosted_global_vector.
reinit(0);
377 ghosted_level_vector.resize(0, 0);
378 solution_ghosted_global_vector.reinit(0);
384template <
typename Number>
389 copy_indices.
clear();
390 solution_copy_indices.clear();
391 copy_indices_global_mine.clear();
392 copy_indices_level_mine.clear();
393 solution_copy_indices_level_mine.clear();
394 component_to_block_map.resize(0);
395 mg_constrained_dofs =
nullptr;
396 ghosted_global_vector.
reinit(0);
397 solution_ghosted_global_vector.reinit(0);
398 ghosted_level_vector.resize(0, 0);
399 perform_plain_copy =
false;
400 perform_renumbered_plain_copy =
false;
401 initialize_dof_vector =
nullptr;
406template <
typename Number>
409 print_indices(std::ostream &os)
const
413 for (
unsigned int i = 0; i < copy_indices[
level].n_cols(); ++i)
414 os <<
"copy_indices[" <<
level <<
"]\t" << copy_indices[
level](0, i)
415 <<
'\t' << copy_indices[
level](1, i) << std::endl;
418 for (
unsigned int level = 0;
level < copy_indices_level_mine.size(); ++
level)
420 for (
unsigned int i = 0; i < copy_indices_level_mine[
level].n_cols(); ++i)
421 os <<
"copy_ifrom [" <<
level <<
"]\t"
422 << copy_indices_level_mine[
level](0, i) <<
'\t'
423 << copy_indices_level_mine[
level](1, i) << std::endl;
425 for (
unsigned int level = 0;
level < copy_indices_global_mine.size(); ++
level)
427 for (
unsigned int i = 0; i < copy_indices_global_mine[
level].n_cols();
429 os <<
"copy_ito [" <<
level <<
"]\t"
430 << copy_indices_global_mine[
level](0, i) <<
'\t'
431 << copy_indices_global_mine[
level](1, i) << std::endl;
437template <
typename Number>
442 std::size_t result =
sizeof(*this);
448 for (
unsigned int i = ghosted_level_vector.min_level();
449 i <= ghosted_level_vector.max_level();
451 result += ghosted_level_vector[i].memory_consumption();
459#include "mg_level_global_transfer.inst"
const IndexSet & locally_owned_mg_dofs(const unsigned int level) const
const Triangulation< dim, spacedim > & get_triangulation() const
const IndexSet & locally_owned_dofs() const
MPI_Comm get_communicator() const
size_type n_elements() const
const std::shared_ptr< const Utilities::MPI::Partitioner > & get_partitioner() const
std::size_t memory_consumption() const
void reinit(const size_type size, const bool omit_zeroing_entries=false)
void print_indices(std::ostream &os) const
std::size_t memory_consumption() const
void fill_and_communicate_copy_indices(const DoFHandler< dim, spacedim > &dof_handler)
virtual unsigned int n_global_levels() const
unsigned int global_to_local(const types::global_dof_index global_index) const
virtual void reinit(const IndexSet &locally_owned_indices, const IndexSet &ghost_indices, const MPI_Comm communicator) override
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define AssertDimension(dim1, dim2)
std::enable_if_t< std::is_fundamental_v< T >, std::size_t > memory_consumption(const T &t)
T max(const T &t, const MPI_Comm mpi_communicator)
T min(const T &t, const MPI_Comm mpi_communicator)
void fill_copy_indices(const DoFHandler< dim, spacedim > &dof_handler, const MGConstrainedDoFs *mg_constrained_dofs, std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > ©_indices, std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > ©_indices_global_mine, std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > ©_indices_level_mine, const bool skip_interface_dofs=true)