37#include <deal.II/multigrid/mg_transfer.templates.h>
48template <
typename VectorType>
49template <
int dim,
int spacedim>
57 copy_indices_global_mine,
58 copy_indices_level_mine);
62 bool my_perform_plain_copy =
68 if (my_perform_plain_copy)
76 for (
const auto ©_index : copy_indices.back())
77 if (copy_index.first != copy_index.second)
79 my_perform_plain_copy =
false;
90 ptria->get_communicator()) == 1);
92 perform_plain_copy = my_perform_plain_copy;
97template <
typename VectorType>
102 copy_indices.clear();
103 copy_indices_global_mine.clear();
104 copy_indices_level_mine.clear();
105 component_to_block_map.resize(0);
106 mg_constrained_dofs =
nullptr;
107 perform_plain_copy =
false;
112template <
typename VectorType>
118 for (
unsigned int i = 0; i < copy_indices[
level].size(); ++i)
119 os <<
"copy_indices[" <<
level <<
"]\t" << copy_indices[
level][i].
first
120 <<
'\t' << copy_indices[
level][i].
second << std::endl;
123 for (
unsigned int level = 0;
level < copy_indices_level_mine.size(); ++
level)
125 for (
unsigned int i = 0; i < copy_indices_level_mine[
level].size(); ++i)
126 os <<
"copy_ifrom [" <<
level <<
"]\t"
127 << copy_indices_level_mine[
level][i].
first <<
'\t'
128 << copy_indices_level_mine[
level][i].
second << std::endl;
130 for (
unsigned int level = 0;
level < copy_indices_global_mine.size(); ++
level)
132 for (
unsigned int i = 0; i < copy_indices_global_mine[
level].size(); ++i)
133 os <<
"copy_ito [" <<
level <<
"]\t"
134 << copy_indices_global_mine[
level][i].
first <<
'\t'
135 << copy_indices_global_mine[
level][i].
second << std::endl;
141template <
typename VectorType>
145 std::size_t result =
sizeof(*this);
160 template <
int dim,
int spacedim,
typename Number>
169 const bool transfer_solution_vectors,
175 &ghosted_level_vector)
179 std::vector<std::pair<types::global_dof_index, types::global_dof_index>>>
182 std::vector<std::pair<types::global_dof_index, types::global_dof_index>>>
183 my_copy_indices_global_mine;
185 std::vector<std::pair<types::global_dof_index, types::global_dof_index>>>
186 my_copy_indices_level_mine;
191 my_copy_indices_global_mine,
192 my_copy_indices_level_mine,
193 !transfer_solution_vectors);
205 std::vector<types::global_dof_index> accessed_indices;
206 ghosted_level_vector.resize(0,
209 std::vector<IndexSet> level_index_set(
214 for (
const auto &indices : my_copy_indices_level_mine[l])
215 accessed_indices.push_back(indices.first);
216 std::vector<types::global_dof_index> accessed_level_indices;
217 for (
const auto &indices : my_copy_indices_global_mine[l])
218 accessed_level_indices.push_back(indices.second);
219 std::sort(accessed_level_indices.begin(), accessed_level_indices.end());
221 level_index_set[l].add_indices(accessed_level_indices.begin(),
222 accessed_level_indices.end());
223 level_index_set[l].compress();
228 std::sort(accessed_indices.begin(), accessed_indices.end());
229 index_set.add_indices(accessed_indices.begin(), accessed_indices.end());
230 index_set.compress();
239 copy_indices_level_mine.resize(
241 copy_indices_global_mine.resize(
243 for (
unsigned int level = 0;
250 *ghosted_level_vector[
level].get_partitioner();
252 auto translate_indices =
253 [&](
const std::vector<
254 std::pair<types::global_dof_index, types::global_dof_index>>
255 & global_copy_indices,
257 local_copy_indices.
reinit(2, global_copy_indices.size());
258 for (
unsigned int i = 0; i < global_copy_indices.size(); ++i)
261 global_copy_indices[i].
first);
263 global_copy_indices[i].
second);
268 translate_indices(my_copy_indices[
level], copy_indices[
level]);
271 translate_indices(my_copy_indices_level_mine[
level],
272 copy_indices_level_mine[
level]);
275 translate_indices(my_copy_indices_global_mine[
level],
276 copy_indices_global_mine[
level]);
281template <
typename Number>
282template <
int dim,
int spacedim>
289 fill_internal(mg_dof,
294 this->copy_indices_global_mine,
295 this->copy_indices_level_mine,
296 ghosted_global_vector,
297 ghosted_level_vector);
303 int have_refinement_edge_dofs = 0;
304 if (mg_constrained_dofs !=
nullptr)
305 for (
unsigned int level = 0;
308 if (mg_constrained_dofs->get_refinement_edge_indices(
level).n_elements() >
311 have_refinement_edge_dofs = 1;
315 fill_internal(mg_dof,
319 this->solution_copy_indices,
320 this->solution_copy_indices_global_mine,
321 this->solution_copy_indices_level_mine,
322 solution_ghosted_global_vector,
323 solution_ghosted_level_vector);
326 this->solution_copy_indices = this->copy_indices;
327 this->solution_copy_indices_global_mine = this->copy_indices_global_mine;
328 this->solution_copy_indices_level_mine = this->copy_indices_level_mine;
329 solution_ghosted_global_vector = ghosted_global_vector;
330 solution_ghosted_level_vector = ghosted_level_vector;
339 const bool my_perform_renumbered_plain_copy =
340 (this->copy_indices.back().n_cols() ==
342 (this->copy_indices_global_mine.back().n_rows() == 0) &&
343 (this->copy_indices_level_mine.back().n_rows() == 0);
345 bool my_perform_plain_copy =
false;
346 if (my_perform_renumbered_plain_copy)
348 my_perform_plain_copy =
true;
352 for (
unsigned int i = 0; i < this->copy_indices.back().n_cols(); ++i)
353 if (this->copy_indices.back()(0, i) != this->copy_indices.back()(1, i))
355 my_perform_plain_copy =
false;
365 perform_renumbered_plain_copy =
370 if (perform_renumbered_plain_copy)
372 for (
unsigned int i = 0; i < this->copy_indices.back().n_cols(); ++i)
375 ghosted_global_vector.
reinit(0);
376 ghosted_level_vector.resize(0, 0);
377 solution_ghosted_global_vector.reinit(0);
378 solution_ghosted_level_vector.resize(0, 0);
384template <
typename Number>
389 copy_indices.
clear();
390 copy_indices_global_mine.clear();
391 copy_indices_level_mine.clear();
392 component_to_block_map.resize(0);
393 mg_constrained_dofs =
nullptr;
394 ghosted_global_vector.
reinit(0);
395 ghosted_level_vector.resize(0, 0);
396 perform_plain_copy =
false;
397 perform_renumbered_plain_copy =
false;
402template <
typename Number>
405 print_indices(std::ostream &os)
const
409 for (
unsigned int i = 0; i < copy_indices[
level].n_cols(); ++i)
410 os <<
"copy_indices[" <<
level <<
"]\t" << copy_indices[
level](0, i)
411 <<
'\t' << copy_indices[
level](1, i) << std::endl;
414 for (
unsigned int level = 0;
level < copy_indices_level_mine.size(); ++
level)
416 for (
unsigned int i = 0; i < copy_indices_level_mine[
level].n_cols(); ++i)
417 os <<
"copy_ifrom [" <<
level <<
"]\t"
418 << copy_indices_level_mine[
level](0, i) <<
'\t'
419 << copy_indices_level_mine[
level](1, i) << std::endl;
421 for (
unsigned int level = 0;
level < copy_indices_global_mine.size(); ++
level)
423 for (
unsigned int i = 0; i < copy_indices_global_mine[
level].n_cols();
425 os <<
"copy_ito [" <<
level <<
"]\t"
426 << copy_indices_global_mine[
level](0, i) <<
'\t'
427 << copy_indices_global_mine[
level](1, i) << std::endl;
433template <
typename Number>
438 std::size_t result =
sizeof(*this);
444 for (
unsigned int i = ghosted_level_vector.min_level();
445 i <= ghosted_level_vector.max_level();
447 result += ghosted_level_vector[i].memory_consumption();
455#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
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 &vector_space_vector_index_set, const IndexSet &read_write_vector_index_set, const MPI_Comm &communicator) override
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define AssertDimension(dim1, dim2)
virtual std::size_t memory_consumption() const override
const std::shared_ptr< const Utilities::MPI::Partitioner > & get_partitioner() const
void reinit(const size_type size, const bool omit_zeroing_entries=false)
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
T min(const T &t, const MPI_Comm &mpi_communicator)
T max(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)