template <int dim, typename Number2, int spacedim>
void
MGLevelGlobalTransfer<LinearAlgebra::distributed::Vector<Number> >::copy_to_mg
-(const DoFHandler<dim,spacedim> &mg_dof_handler,
+(const DoFHandler<dim,spacedim> &mg_dof_handler,
MGLevelObject<LinearAlgebra::distributed::Vector<Number> > &dst,
const LinearAlgebra::distributed::Vector<Number2> &src,
const bool solution_transfer) const
if (perform_plain_copy)
{
- // In this case, we can simply copy the local range (in parallel by
- // VectorView)
+ // In this case, we can simply copy the local range.
AssertDimension(dst[dst.max_level()].local_size(), src.local_size());
- VectorView<Number> dst_view (src.local_size(), dst[dst.max_level()].begin());
- VectorView<Number2> src_view (src.local_size(), src.begin());
- static_cast<Vector<Number> &>(dst_view) = static_cast<Vector<Number2> &>(src_view);
+ dst[dst.max_level()].copy_locally_owned_data_from(src);
+ return;
+ }
+ else if (perform_renumbered_plain_copy)
+ {
+ Assert(copy_indices_level_mine.back().empty(), ExcInternalError());
+ LinearAlgebra::distributed::Vector<Number> &dst_level = dst[dst.max_level()];
+ for (std::pair<unsigned int,unsigned int> i : this_copy_indices.back())
+ dst_level.local_element(i.second) = src.local_element(i.first);
return;
}
template <int dim, typename Number2, int spacedim>
void
MGLevelGlobalTransfer<LinearAlgebra::distributed::Vector<Number> >::copy_from_mg
-(const DoFHandler<dim,spacedim> &mg_dof_handler,
+(const DoFHandler<dim,spacedim> &mg_dof_handler,
LinearAlgebra::distributed::Vector<Number2> &dst,
const MGLevelObject<LinearAlgebra::distributed::Vector<Number> > &src) const
{
AssertIndexRange(src.min_level(), src.max_level()+1);
if (perform_plain_copy)
{
- // In this case, we can simply copy the local range (in parallel by
- // VectorView). To avoid having stray data in ghost entries of the
- // destination, make sure to clear them here.
+ // In this case, we can simply copy the local range. To avoid having
+ // stray data in ghost entries of the destination, make sure to clear
+ // them here.
+ dst.zero_out_ghosts();
+ AssertDimension(src[src.max_level()].local_size(), dst.local_size());
+ dst.copy_locally_owned_data_from(src[src.max_level()]);
+ return;
+ }
+ else if (perform_renumbered_plain_copy)
+ {
+ Assert(copy_indices_global_mine.back().empty(), ExcInternalError());
+ const LinearAlgebra::distributed::Vector<Number> &src_level = src[src.max_level()];
dst.zero_out_ghosts();
- AssertDimension(dst.local_size(), src[src.max_level()].local_size());
- VectorView<Number2> dst_view (dst.local_size(), dst.begin());
- VectorView<Number> src_view (dst.local_size(), src[src.max_level()].begin());
- static_cast<Vector<Number2> &>(dst_view) = static_cast<Vector<Number> &>(src_view);
+ for (std::pair<unsigned int,unsigned int> i : copy_indices.back())
+ dst.local_element(i.first) = src_level.local_element(i.second);
return;
}
template <int dim, typename Number2, int spacedim>
void
MGLevelGlobalTransfer<LinearAlgebra::distributed::Vector<Number> >::copy_from_mg_add
-(const DoFHandler<dim,spacedim> &/*mg_dof_handler*/,
+(const DoFHandler<dim,spacedim> &/*mg_dof_handler*/,
LinearAlgebra::distributed::Vector<Number2> &dst,
const MGLevelObject<LinearAlgebra::distributed::Vector<Number> > &src) const
{
copy_indices_global_mine, copy_indices_level_mine);
// check if we can run a plain copy operation between the global DoFs and
- // the finest level, on the current processor.
+ // the finest level.
bool
my_perform_plain_copy =
(copy_indices.back().size() == mg_dof.locally_owned_dofs().n_elements())
copy_indices_level_mine.clear();
component_to_block_map.resize(0);
mg_constrained_dofs = nullptr;
+ perform_plain_copy = false;
}
solution_ghosted_global_vector,
solution_ghosted_level_vector);
- perform_plain_copy = this->copy_indices.back().size()
- == mg_dof.locally_owned_dofs().n_elements();
- if (perform_plain_copy)
+ bool my_perform_renumbered_plain_copy = (this->copy_indices.back().size() ==
+ mg_dof.locally_owned_dofs().n_elements());
+ bool my_perform_plain_copy = false;
+ if (my_perform_renumbered_plain_copy)
{
+ my_perform_plain_copy = true;
AssertDimension(this->copy_indices_global_mine.back().size(), 0);
AssertDimension(this->copy_indices_level_mine.back().size(), 0);
if (this->copy_indices.back()[i].first !=
this->copy_indices.back()[i].second)
{
- perform_plain_copy = false;
+ my_perform_plain_copy = false;
break;
}
}
+
+ // now do a global reduction over all processors to see what operation
+ // they can agree upon
perform_plain_copy =
- Utilities::MPI::min(static_cast<int>(perform_plain_copy),
+ Utilities::MPI::min(static_cast<int>(my_perform_plain_copy),
+ mpi_communicator);
+ perform_renumbered_plain_copy =
+ Utilities::MPI::min(static_cast<int>(my_perform_renumbered_plain_copy),
mpi_communicator);
// if we do a plain copy, no need to hold additional ghosted vectors
- if (perform_plain_copy)
+ if (perform_renumbered_plain_copy)
{
ghosted_global_vector.reinit(0);
ghosted_level_vector.resize(0, 0);
mg_constrained_dofs = nullptr;
ghosted_global_vector.reinit(0);
ghosted_level_vector.resize(0, 0);
+ perform_plain_copy = false;
+ perform_renumbered_plain_copy = false;
}