From 29de9ecfa02b230a9a826ba355b60fd3a118026b Mon Sep 17 00:00:00 2001 From: Peter Munch Date: Sun, 24 Apr 2022 22:46:03 +0200 Subject: [PATCH] Introduce MGTransferBlockGlobalCoarsening --- include/deal.II/base/mg_level_object.h | 14 ++ .../multigrid/mg_transfer_global_coarsening.h | 33 +++ .../mg_transfer_global_coarsening.templates.h | 27 +++ .../multigrid/mg_transfer_matrix_free.h | 39 +++- .../mg_transfer_matrix_free.templates.h | 86 ++++++++ .../mg_transfer_global_coarsening.inst.in | 8 + source/multigrid/mg_transfer_matrix_free.cc | 86 ++------ tests/multigrid-global-coarsening/block_01.cc | 190 ++++++++++++++++++ .../block_01.with_p4est=true.mpirun=1.output | 73 +++++++ 9 files changed, 473 insertions(+), 83 deletions(-) create mode 100644 include/deal.II/multigrid/mg_transfer_matrix_free.templates.h create mode 100644 tests/multigrid-global-coarsening/block_01.cc create mode 100644 tests/multigrid-global-coarsening/block_01.with_p4est=true.mpirun=1.output diff --git a/include/deal.II/base/mg_level_object.h b/include/deal.II/base/mg_level_object.h index 88a9e3ae47..fd10610f5d 100644 --- a/include/deal.II/base/mg_level_object.h +++ b/include/deal.II/base/mg_level_object.h @@ -97,6 +97,12 @@ public: const Object & operator[](const unsigned int level) const; + /** + * Return object on level max. + */ + const Object & + back() const; + /** * Delete all previous contents of this object and reset its size according * to the values of @p new_minlevel and @p new_maxlevel. @@ -230,6 +236,14 @@ MGLevelObject::operator[](const unsigned int i) const } +template +const Object & +MGLevelObject::back() const +{ + return this->operator[](this->max_level()); +} + + template template void diff --git a/include/deal.II/multigrid/mg_transfer_global_coarsening.h b/include/deal.II/multigrid/mg_transfer_global_coarsening.h index 7b93d0e367..8948f5ba39 100644 --- a/include/deal.II/multigrid/mg_transfer_global_coarsening.h +++ b/include/deal.II/multigrid/mg_transfer_global_coarsening.h @@ -28,6 +28,7 @@ #include #include +#include DEAL_II_NAMESPACE_OPEN @@ -616,6 +617,38 @@ private: +/** + * This class works with LinearAlgebra::distributed::BlockVector and + * performs exactly the same transfer operations for each block as + * MGTransferGlobalCoarsening. + */ +template +class MGTransferBlockGlobalCoarsening + : public MGTransferBlockMatrixFreeBase< + dim, + typename VectorType::value_type, + MGTransferGlobalCoarsening> +{ +public: + /** + * Constructor. + */ + MGTransferBlockGlobalCoarsening( + const MGTransferGlobalCoarsening &transfer_operator); + +protected: + const MGTransferGlobalCoarsening & + get_matrix_free_transfer(const unsigned int b) const override; + +private: + /** + * Non-block version of transfer operation. + */ + const MGTransferGlobalCoarsening &transfer_operator; +}; + + + #ifndef DOXYGEN /* ----------------------- Inline functions --------------------------------- */ diff --git a/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h b/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h index 3a14658609..70ab1f48b4 100644 --- a/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h +++ b/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h @@ -45,6 +45,7 @@ #include #include +#include DEAL_II_NAMESPACE_OPEN @@ -3163,6 +3164,32 @@ MGTwoLevelTransfer>:: return size; } + + +template +MGTransferBlockGlobalCoarsening:: + MGTransferBlockGlobalCoarsening( + const MGTransferGlobalCoarsening &transfer_operator) + : MGTransferBlockMatrixFreeBase>( + true) + , transfer_operator(transfer_operator) +{} + + + +template +const MGTransferGlobalCoarsening & +MGTransferBlockGlobalCoarsening::get_matrix_free_transfer( + const unsigned int b) const +{ + (void)b; + AssertDimension(b, 0); + return transfer_operator; +} + + DEAL_II_NAMESPACE_CLOSE #endif diff --git a/include/deal.II/multigrid/mg_transfer_matrix_free.h b/include/deal.II/multigrid/mg_transfer_matrix_free.h index 3368326066..df6a6b9cc1 100644 --- a/include/deal.II/multigrid/mg_transfer_matrix_free.h +++ b/include/deal.II/multigrid/mg_transfer_matrix_free.h @@ -321,6 +321,10 @@ class MGTransferBlockMatrixFreeBase : public MGTransferBase> { public: + MGTransferBlockMatrixFreeBase(const bool same_for_all) + : same_for_all(same_for_all) + {} + /** * Prolongate a vector from level to_level-1 to level * to_level using the embedding matrices of the underlying finite @@ -428,15 +432,17 @@ public: protected: /** - * Non-block matrix-free versions of transfer operation. + * Return the right non-block transfer operator. Has to be implemented by + * the derived class. */ - std::vector matrix_free_transfer_vector; + virtual const TransferType & + get_matrix_free_transfer(const unsigned int b) const = 0; /** * A flag to indicate whether the same DoFHandler is used for all * the components or if each block has its own DoFHandler. */ - bool same_for_all; + const bool same_for_all; }; @@ -520,6 +526,16 @@ public: */ std::size_t memory_consumption() const; + +protected: + const MGTransferMatrixFree & + get_matrix_free_transfer(const unsigned int b) const override; + +private: + /** + * Non-block matrix-free versions of transfer operation. + */ + std::vector> matrix_free_transfer_vector; }; @@ -631,7 +647,6 @@ MGTransferBlockMatrixFreeBase::copy_to_mg( MGLevelObject> &dst, const LinearAlgebra::distributed::BlockVector & src) const { - AssertDimension(matrix_free_transfer_vector.size(), 1); Assert(same_for_all, ExcMessage( "This object was initialized with support for usage with one " @@ -673,9 +688,8 @@ MGTransferBlockMatrixFreeBase::copy_to_mg( for (unsigned int b = 0; b < n_blocks; ++b) { const unsigned int data_block = same_for_all ? 0 : b; - matrix_free_transfer_vector[data_block].copy_to_mg(*dof_handler[b], - dst_non_block, - src.block(b)); + get_matrix_free_transfer(data_block) + .copy_to_mg(*dof_handler[b], dst_non_block, src.block(b)); for (unsigned int l = min_level; l <= max_level; ++l) dst[l].block(b) = dst_non_block[l]; @@ -694,7 +708,11 @@ MGTransferBlockMatrixFreeBase::copy_from_mg( const MGLevelObject> &src) const { - AssertDimension(matrix_free_transfer_vector.size(), 1); + Assert(same_for_all, + ExcMessage( + "This object was initialized with support for usage with one " + "DoFHandler for each block, but this method assumes that " + "the same DoFHandler is used for all the blocks!")); const std::vector *> mg_dofs(dst.n_blocks(), &dof_handler); @@ -734,9 +752,8 @@ MGTransferBlockMatrixFreeBase::copy_from_mg( src_non_block[l] = src[l].block(b); } const unsigned int data_block = same_for_all ? 0 : b; - matrix_free_transfer_vector[data_block].copy_from_mg(*dof_handler[b], - dst.block(b), - src_non_block); + get_matrix_free_transfer(data_block) + .copy_from_mg(*dof_handler[b], dst.block(b), src_non_block); } } diff --git a/include/deal.II/multigrid/mg_transfer_matrix_free.templates.h b/include/deal.II/multigrid/mg_transfer_matrix_free.templates.h new file mode 100644 index 0000000000..e45d295d92 --- /dev/null +++ b/include/deal.II/multigrid/mg_transfer_matrix_free.templates.h @@ -0,0 +1,86 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2016 - 2021 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE.md at +// the top level directory of deal.II. +// +// --------------------------------------------------------------------- + +#ifndef dealii_mg_transfer_matrix_free_templates_h +#define dealii_mg_transfer_matrix_free_templates_h + +#include + +#include + + +DEAL_II_NAMESPACE_OPEN + +template +void +MGTransferBlockMatrixFreeBase::prolongate( + const unsigned int to_level, + LinearAlgebra::distributed::BlockVector & dst, + const LinearAlgebra::distributed::BlockVector &src) const +{ + const unsigned int n_blocks = src.n_blocks(); + AssertDimension(dst.n_blocks(), n_blocks); + + for (unsigned int b = 0; b < n_blocks; ++b) + { + const unsigned int data_block = same_for_all ? 0 : b; + get_matrix_free_transfer(data_block) + .prolongate(to_level, dst.block(b), src.block(b)); + } +} + + + +template +void +MGTransferBlockMatrixFreeBase::prolongate_and_add( + const unsigned int to_level, + LinearAlgebra::distributed::BlockVector & dst, + const LinearAlgebra::distributed::BlockVector &src) const +{ + const unsigned int n_blocks = src.n_blocks(); + AssertDimension(dst.n_blocks(), n_blocks); + + for (unsigned int b = 0; b < n_blocks; ++b) + { + const unsigned int data_block = same_for_all ? 0 : b; + get_matrix_free_transfer(data_block) + .prolongate_and_add(to_level, dst.block(b), src.block(b)); + } +} + + + +template +void +MGTransferBlockMatrixFreeBase::restrict_and_add( + const unsigned int from_level, + LinearAlgebra::distributed::BlockVector & dst, + const LinearAlgebra::distributed::BlockVector &src) const +{ + const unsigned int n_blocks = src.n_blocks(); + AssertDimension(dst.n_blocks(), n_blocks); + + for (unsigned int b = 0; b < n_blocks; ++b) + { + const unsigned int data_block = same_for_all ? 0 : b; + get_matrix_free_transfer(data_block) + .restrict_and_add(from_level, dst.block(b), src.block(b)); + } +} + +DEAL_II_NAMESPACE_CLOSE + +#endif diff --git a/source/multigrid/mg_transfer_global_coarsening.inst.in b/source/multigrid/mg_transfer_global_coarsening.inst.in index e9a9eb71d7..209a7ee8c1 100644 --- a/source/multigrid/mg_transfer_global_coarsening.inst.in +++ b/source/multigrid/mg_transfer_global_coarsening.inst.in @@ -19,6 +19,14 @@ for (deal_II_dimension : DIMENSIONS; S1 : REAL_SCALARS) { template class MGTwoLevelTransfer>; + template class MGTransferBlockMatrixFreeBase< + deal_II_dimension, + S1, + MGTransferGlobalCoarsening>>; + template class MGTransferBlockGlobalCoarsening< + deal_II_dimension, + LinearAlgebra::distributed::Vector>; } for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS) diff --git a/source/multigrid/mg_transfer_matrix_free.cc b/source/multigrid/mg_transfer_matrix_free.cc index 19fcb4295f..af000da93f 100644 --- a/source/multigrid/mg_transfer_matrix_free.cc +++ b/source/multigrid/mg_transfer_matrix_free.cc @@ -33,7 +33,7 @@ #include #include -#include +#include #include @@ -707,8 +707,10 @@ MGTransferMatrixFree::memory_consumption() const template MGTransferBlockMatrixFree::MGTransferBlockMatrixFree( const MGConstrainedDoFs &mg_c) + : MGTransferBlockMatrixFreeBase>(true) { - this->same_for_all = true; this->matrix_free_transfer_vector.emplace_back(mg_c); } @@ -717,8 +719,10 @@ MGTransferBlockMatrixFree::MGTransferBlockMatrixFree( template MGTransferBlockMatrixFree::MGTransferBlockMatrixFree( const std::vector &mg_c) + : MGTransferBlockMatrixFreeBase>(false) { - this->same_for_all = false; for (const auto &constrained_block_dofs : mg_c) this->matrix_free_transfer_vector.emplace_back(constrained_block_dofs); } @@ -795,82 +799,20 @@ MGTransferBlockMatrixFree::memory_consumption() const template -void -MGTransferBlockMatrixFree::clear() -{ - this->matrix_free_transfer_vector.clear(); -} - - - -template -void -MGTransferBlockMatrixFreeBase::prolongate( - const unsigned int to_level, - LinearAlgebra::distributed::BlockVector & dst, - const LinearAlgebra::distributed::BlockVector &src) const -{ - const unsigned int n_blocks = src.n_blocks(); - AssertDimension(dst.n_blocks(), n_blocks); - - if (!same_for_all) - AssertDimension(matrix_free_transfer_vector.size(), n_blocks); - - for (unsigned int b = 0; b < n_blocks; ++b) - { - const unsigned int data_block = same_for_all ? 0 : b; - matrix_free_transfer_vector[data_block].prolongate(to_level, - dst.block(b), - src.block(b)); - } -} - - - -template -void -MGTransferBlockMatrixFreeBase::prolongate_and_add( - const unsigned int to_level, - LinearAlgebra::distributed::BlockVector & dst, - const LinearAlgebra::distributed::BlockVector &src) const +const MGTransferMatrixFree & +MGTransferBlockMatrixFree::get_matrix_free_transfer( + const unsigned int b) const { - const unsigned int n_blocks = src.n_blocks(); - AssertDimension(dst.n_blocks(), n_blocks); - - if (!same_for_all) - AssertDimension(matrix_free_transfer_vector.size(), n_blocks); - - for (unsigned int b = 0; b < n_blocks; ++b) - { - const unsigned int data_block = same_for_all ? 0 : b; - matrix_free_transfer_vector[data_block].prolongate_and_add(to_level, - dst.block(b), - src.block(b)); - } + return matrix_free_transfer_vector[b]; } -template +template void -MGTransferBlockMatrixFreeBase::restrict_and_add( - const unsigned int from_level, - LinearAlgebra::distributed::BlockVector & dst, - const LinearAlgebra::distributed::BlockVector &src) const +MGTransferBlockMatrixFree::clear() { - const unsigned int n_blocks = src.n_blocks(); - AssertDimension(dst.n_blocks(), n_blocks); - - if (!same_for_all) - AssertDimension(matrix_free_transfer_vector.size(), n_blocks); - - for (unsigned int b = 0; b < n_blocks; ++b) - { - const unsigned int data_block = same_for_all ? 0 : b; - matrix_free_transfer_vector[data_block].restrict_and_add(from_level, - dst.block(b), - src.block(b)); - } + this->matrix_free_transfer_vector.clear(); } diff --git a/tests/multigrid-global-coarsening/block_01.cc b/tests/multigrid-global-coarsening/block_01.cc new file mode 100644 index 0000000000..32316620a4 --- /dev/null +++ b/tests/multigrid-global-coarsening/block_01.cc @@ -0,0 +1,190 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2022 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE.md at +// the top level directory of deal.II. +// +// --------------------------------------------------------------------- + + +/** + * Test global-coarsening multigrid for block vectors. + */ + +#include +#include +#include +#include + +#include + +#include +#include + +#include +#include +#include +#include +#include + +#include +#include + +#include +#include + +#include +#include + +#include "../tests.h" + +template +void +test(const unsigned int n_refinements, const unsigned int fe_degree) +{ + using VectorType = LinearAlgebra::distributed::Vector; + + const unsigned int min_level = 0; + const unsigned int max_level = n_refinements; + + parallel::distributed::Triangulation tria(MPI_COMM_WORLD); + + GridGenerator::hyper_cube(tria); + tria.refine_global(n_refinements); + + const auto trias = + MGTransferGlobalCoarseningTools::create_geometric_coarsening_sequence(tria); + + MGLevelObject> dof_handlers(min_level, max_level); + MGLevelObject> dof_handlers_system(min_level, max_level); + MGLevelObject> + partitioners(min_level, max_level); + MGLevelObject> + partitioners_system(min_level, max_level); + MGLevelObject> transfers(min_level, + max_level); + MGLevelObject> transfers_system( + min_level, max_level); + + for (auto l = min_level; l <= max_level; ++l) + { + dof_handlers[l].reinit(*trias[l]); + dof_handlers[l].distribute_dofs(FE_Q(fe_degree)); + partitioners[l] = std::make_shared( + dof_handlers[l].locally_owned_dofs(), + DoFTools::extract_locally_relevant_dofs(dof_handlers[l]), + MPI_COMM_WORLD); + + + dof_handlers_system[l].reinit(*trias[l]); + dof_handlers_system[l].distribute_dofs( + FESystem(FE_Q(fe_degree), dim)); + DoFRenumbering::block_wise(dof_handlers_system[l]); + partitioners_system[l] = + std::make_shared( + dof_handlers_system[l].locally_owned_dofs(), + DoFTools::extract_locally_relevant_dofs(dof_handlers_system[l]), + MPI_COMM_WORLD); + } + + for (unsigned int l = min_level; l < max_level; ++l) + { + transfers[l + 1].reinit(dof_handlers[l + 1], dof_handlers[l]); + transfers_system[l + 1].reinit(dof_handlers_system[l + 1], + dof_handlers_system[l]); + } + + + + MGTransferGlobalCoarsening transfer( + transfers, [&](const auto l, auto &vec) { vec.reinit(partitioners[l]); }); + + MGTransferGlobalCoarsening transfer_system( + transfers_system, + [&](const auto l, auto &vec) { vec.reinit(partitioners_system[l]); }); + + MGTransferBlockGlobalCoarsening transfer_block(transfer); + + LinearAlgebra::distributed::BlockVector vec_block(dim); + for (unsigned int d = 0; d < dim; ++d) + vec_block.block(d).reinit(partitioners.back()); + MGLevelObject> + vec_levels_block(min_level, max_level); + + LinearAlgebra::distributed::Vector vec_system( + partitioners_system.back()); + MGLevelObject> vec_levels_system( + min_level, max_level); + + for (unsigned int b = 0, c = 0; b < dim; ++b) + for (unsigned int i = 0; i < dof_handlers.back().n_dofs(); ++i, ++c) + { + vec_block.block(b)[i] = c; + vec_system[c] = c; + } + + // test copy_to_mg() + transfer_block.copy_to_mg(dof_handlers.back(), vec_levels_block, vec_block); + transfer_system.copy_to_mg(dof_handlers_system.back(), + vec_levels_system, + vec_system); + + vec_levels_block.back().print(deallog.get_file_stream()); + vec_levels_system.back().print(deallog.get_file_stream()); + + // test copy_from_mg() + transfer_block.copy_from_mg(dof_handlers.back(), vec_block, vec_levels_block); + transfer_system.copy_from_mg(dof_handlers_system.back(), + vec_system, + vec_levels_system); + + vec_block.print(deallog.get_file_stream()); + vec_system.print(deallog.get_file_stream()); + + + // test restrict_and_add() + for (unsigned int l = max_level; l > min_level; --l) + { + transfer_block.restrict_and_add(l, + vec_levels_block[l - 1], + vec_levels_block[l]); + transfer_system.restrict_and_add(l, + vec_levels_system[l - 1], + vec_levels_system[l]); + + vec_levels_block[l - 1].print(deallog.get_file_stream()); + vec_levels_system[l - 1].print(deallog.get_file_stream()); + } + + // test prolongate() + for (unsigned int l = min_level; l < max_level; ++l) + { + transfer_block.prolongate(l + 1, + vec_levels_block[l + 1], + vec_levels_block[l]); + transfer_system.prolongate(l + 1, + vec_levels_system[l + 1], + vec_levels_system[l]); + + vec_levels_block[l + 1].print(deallog.get_file_stream()); + vec_levels_system[l + 1].print(deallog.get_file_stream()); + } +} + +int +main(int argc, char **argv) +{ + Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1); + MPILogInitAll all; + + deallog.precision(8); + + test<2>(2, 1); +} diff --git a/tests/multigrid-global-coarsening/block_01.with_p4est=true.mpirun=1.output b/tests/multigrid-global-coarsening/block_01.with_p4est=true.mpirun=1.output new file mode 100644 index 0000000000..29dae0a482 --- /dev/null +++ b/tests/multigrid-global-coarsening/block_01.with_p4est=true.mpirun=1.output @@ -0,0 +1,73 @@ + +Process #0 +Local range: [0, 25), global size: 25 +Vector data: +0.000e+00 1.000e+00 2.000e+00 3.000e+00 4.000e+00 5.000e+00 6.000e+00 7.000e+00 8.000e+00 9.000e+00 1.000e+01 1.100e+01 1.200e+01 1.300e+01 1.400e+01 1.500e+01 1.600e+01 1.700e+01 1.800e+01 1.900e+01 2.000e+01 2.100e+01 2.200e+01 2.300e+01 2.400e+01 +Process #0 +Local range: [0, 25), global size: 25 +Vector data: +2.500e+01 2.600e+01 2.700e+01 2.800e+01 2.900e+01 3.000e+01 3.100e+01 3.200e+01 3.300e+01 3.400e+01 3.500e+01 3.600e+01 3.700e+01 3.800e+01 3.900e+01 4.000e+01 4.100e+01 4.200e+01 4.300e+01 4.400e+01 4.500e+01 4.600e+01 4.700e+01 4.800e+01 4.900e+01 +Process #0 +Local range: [0, 50), global size: 50 +Vector data: +0.000e+00 1.000e+00 2.000e+00 3.000e+00 4.000e+00 5.000e+00 6.000e+00 7.000e+00 8.000e+00 9.000e+00 1.000e+01 1.100e+01 1.200e+01 1.300e+01 1.400e+01 1.500e+01 1.600e+01 1.700e+01 1.800e+01 1.900e+01 2.000e+01 2.100e+01 2.200e+01 2.300e+01 2.400e+01 2.500e+01 2.600e+01 2.700e+01 2.800e+01 2.900e+01 3.000e+01 3.100e+01 3.200e+01 3.300e+01 3.400e+01 3.500e+01 3.600e+01 3.700e+01 3.800e+01 3.900e+01 4.000e+01 4.100e+01 4.200e+01 4.300e+01 4.400e+01 4.500e+01 4.600e+01 4.700e+01 4.800e+01 4.900e+01 +Process #0 +Local range: [0, 25), global size: 25 +Vector data: +0.000e+00 1.000e+00 2.000e+00 3.000e+00 4.000e+00 5.000e+00 6.000e+00 7.000e+00 8.000e+00 9.000e+00 1.000e+01 1.100e+01 1.200e+01 1.300e+01 1.400e+01 1.500e+01 1.600e+01 1.700e+01 1.800e+01 1.900e+01 2.000e+01 2.100e+01 2.200e+01 2.300e+01 2.400e+01 +Process #0 +Local range: [0, 25), global size: 25 +Vector data: +2.500e+01 2.600e+01 2.700e+01 2.800e+01 2.900e+01 3.000e+01 3.100e+01 3.200e+01 3.300e+01 3.400e+01 3.500e+01 3.600e+01 3.700e+01 3.800e+01 3.900e+01 4.000e+01 4.100e+01 4.200e+01 4.300e+01 4.400e+01 4.500e+01 4.600e+01 4.700e+01 4.800e+01 4.900e+01 +Process #0 +Local range: [0, 50), global size: 50 +Vector data: +0.000e+00 1.000e+00 2.000e+00 3.000e+00 4.000e+00 5.000e+00 6.000e+00 7.000e+00 8.000e+00 9.000e+00 1.000e+01 1.100e+01 1.200e+01 1.300e+01 1.400e+01 1.500e+01 1.600e+01 1.700e+01 1.800e+01 1.900e+01 2.000e+01 2.100e+01 2.200e+01 2.300e+01 2.400e+01 2.500e+01 2.600e+01 2.700e+01 2.800e+01 2.900e+01 3.000e+01 3.100e+01 3.200e+01 3.300e+01 3.400e+01 3.500e+01 3.600e+01 3.700e+01 3.800e+01 3.900e+01 4.000e+01 4.100e+01 4.200e+01 4.300e+01 4.400e+01 4.500e+01 4.600e+01 4.700e+01 4.800e+01 4.900e+01 +Process #0 +Local range: [0, 9), global size: 9 +Vector data: +2.250e+00 1.475e+01 2.275e+01 4.150e+01 2.400e+01 4.525e+01 3.900e+01 5.875e+01 5.175e+01 +Process #0 +Local range: [0, 9), global size: 9 +Vector data: +5.850e+01 8.975e+01 9.775e+01 1.415e+02 8.025e+01 1.202e+02 9.525e+01 1.338e+02 1.080e+02 +Process #0 +Local range: [0, 18), global size: 18 +Vector data: +2.250e+00 1.475e+01 2.275e+01 4.150e+01 2.400e+01 4.525e+01 3.900e+01 5.875e+01 5.175e+01 5.850e+01 8.975e+01 9.775e+01 1.415e+02 8.025e+01 1.202e+02 9.525e+01 1.338e+02 1.080e+02 +Process #0 +Local range: [0, 4), global size: 4 +Vector data: +3.138e+01 6.438e+01 9.012e+01 1.141e+02 +Process #0 +Local range: [0, 4), global size: 4 +Vector data: +1.876e+02 2.206e+02 2.464e+02 2.704e+02 +Process #0 +Local range: [0, 8), global size: 8 +Vector data: +3.138e+01 6.438e+01 9.012e+01 1.141e+02 1.876e+02 2.206e+02 2.464e+02 2.704e+02 +Process #0 +Local range: [0, 9), global size: 9 +Vector data: +3.138e+01 4.788e+01 6.075e+01 7.500e+01 6.438e+01 8.925e+01 9.012e+01 1.021e+02 1.141e+02 +Process #0 +Local range: [0, 9), global size: 9 +Vector data: +1.876e+02 2.041e+02 2.170e+02 2.312e+02 2.206e+02 2.455e+02 2.464e+02 2.584e+02 2.704e+02 +Process #0 +Local range: [0, 18), global size: 18 +Vector data: +3.138e+01 4.788e+01 6.075e+01 7.500e+01 6.438e+01 8.925e+01 9.012e+01 1.021e+02 1.141e+02 1.876e+02 2.041e+02 2.170e+02 2.312e+02 2.206e+02 2.455e+02 2.464e+02 2.584e+02 2.704e+02 +Process #0 +Local range: [0, 25), global size: 25 +Vector data: +3.138e+01 3.962e+01 4.606e+01 5.375e+01 4.788e+01 6.144e+01 6.075e+01 6.788e+01 7.500e+01 5.612e+01 6.912e+01 6.438e+01 7.681e+01 8.212e+01 8.925e+01 7.544e+01 8.200e+01 8.856e+01 9.012e+01 9.612e+01 1.021e+02 9.512e+01 1.017e+02 1.081e+02 1.141e+02 +Process #0 +Local range: [0, 25), global size: 25 +Vector data: +1.876e+02 1.959e+02 2.023e+02 2.100e+02 2.041e+02 2.177e+02 2.170e+02 2.241e+02 2.312e+02 2.124e+02 2.254e+02 2.206e+02 2.331e+02 2.384e+02 2.455e+02 2.317e+02 2.382e+02 2.448e+02 2.464e+02 2.524e+02 2.584e+02 2.514e+02 2.579e+02 2.644e+02 2.704e+02 +Process #0 +Local range: [0, 50), global size: 50 +Vector data: +3.138e+01 3.962e+01 4.606e+01 5.375e+01 4.788e+01 6.144e+01 6.075e+01 6.788e+01 7.500e+01 5.612e+01 6.912e+01 6.438e+01 7.681e+01 8.212e+01 8.925e+01 7.544e+01 8.200e+01 8.856e+01 9.012e+01 9.612e+01 1.021e+02 9.512e+01 1.017e+02 1.081e+02 1.141e+02 1.876e+02 1.959e+02 2.023e+02 2.100e+02 2.041e+02 2.177e+02 2.170e+02 2.241e+02 2.312e+02 2.124e+02 2.254e+02 2.206e+02 2.331e+02 2.384e+02 2.455e+02 2.317e+02 2.382e+02 2.448e+02 2.464e+02 2.524e+02 2.584e+02 2.514e+02 2.579e+02 2.644e+02 2.704e+02 -- 2.39.5