From 01b94295b61e1866a386238beeb7f38018ac5420 Mon Sep 17 00:00:00 2001 From: David Wells Date: Wed, 4 Oct 2023 11:28:55 -0400 Subject: [PATCH] Add code to serialize PETSc vectors. --- include/deal.II/lac/petsc_vector_base.h | 96 +++++++++++++++++++ tests/petsc/serialize_vector_01.cc | 63 ++++++++++++ .../petsc/serialize_vector_01.mpirun=2.output | 9 ++ 3 files changed, 168 insertions(+) create mode 100644 tests/petsc/serialize_vector_01.cc create mode 100644 tests/petsc/serialize_vector_01.mpirun=2.output diff --git a/include/deal.II/lac/petsc_vector_base.h b/include/deal.II/lac/petsc_vector_base.h index a2dec0f63e..d936d92aa5 100644 --- a/include/deal.II/lac/petsc_vector_base.h +++ b/include/deal.II/lac/petsc_vector_base.h @@ -28,6 +28,8 @@ # include # include +# include + # include # include @@ -724,6 +726,41 @@ namespace PETScWrappers const bool scientific = true, const bool across = true) const; + /** + * Write the data of this object to a stream for the purpose of + * serialization using the [BOOST serialization + * library](https://www.boost.org/doc/libs/1_74_0/libs/serialization/doc/index.html). + * + * @note Each processor only serializes its own locally owned values. + */ + template + void + save(Archive &ar, const unsigned int version) const; + + /** + * Read the data of this object from a stream for the purpose of + * serialization using the [BOOST serialization + * library](https://www.boost.org/doc/libs/1_74_0/libs/serialization/doc/index.html). + */ + template + void + load(Archive &ar, const unsigned int version); + +# ifdef DOXYGEN + /** + * Write and read the data of this object from a stream for the purpose + * of serialization using the [BOOST serialization + * library](https://www.boost.org/doc/libs/1_74_0/libs/serialization/doc/index.html). + */ + template + void + serialize(Archive &archive, const unsigned int version); +# else + // This macro defines the serialize() method that is compatible with + // the templated save() and load() method that have been implemented. + BOOST_SERIALIZATION_SPLIT_MEMBER() +# endif + /** * Swap the contents of this vector and the other vector @p v. One could * do this operation with a temporary variable and copying over the data @@ -1271,6 +1308,65 @@ namespace PETScWrappers } } + template + inline void + VectorBase::save(Archive &ar, const unsigned int) const + { + // forward to serialization function in the base class. + ar &static_cast(*this); + ar &size(); + ar &local_range(); + + const PetscScalar *array = nullptr; + int ierr = VecGetArrayRead(*this, &array); + AssertThrow(ierr == 0, ExcPETScError(ierr)); + + boost::serialization::array_wrapper wrapper( + array, locally_owned_size()); + ar &wrapper; + + ierr = VecRestoreArrayRead(*this, &array); + AssertThrow(ierr == 0, ExcPETScError(ierr)); + } + + + + template + inline void + VectorBase::load(Archive &ar, const unsigned int) + { + ar &static_cast(*this); + + size_type size = 0; + std::pair local_range; + + ar &size; + Assert(size == this->size(), + ExcMessage("The serialized value of size (" + std::to_string(size) + + ") does not match the current size (" + + std::to_string(this->size()) + ")")); + (void)size; + ar &local_range; + Assert(local_range == this->local_range(), + ExcMessage("The serialized value of local_range (" + + std::to_string(local_range.first) + ", " + + std::to_string(local_range.second) + + ") does not match the current local_range (" + + std::to_string(this->local_range().first) + ", " + + std::to_string(this->local_range().second) + ")")); + (void)local_range; + + PetscScalar *array = nullptr; + int ierr = VecGetArray(petsc_vector(), &array); + AssertThrow(ierr == 0, ExcPETScError(ierr)); + + boost::serialization::array_wrapper velocity_wrapper( + array, locally_owned_size()); + ar &velocity_wrapper; + + ierr = VecRestoreArray(petsc_vector(), &array); + AssertThrow(ierr == 0, ExcPETScError(ierr)); + } # endif // DOXYGEN } // namespace PETScWrappers diff --git a/tests/petsc/serialize_vector_01.cc b/tests/petsc/serialize_vector_01.cc new file mode 100644 index 0000000000..5a51af471c --- /dev/null +++ b/tests/petsc/serialize_vector_01.cc @@ -0,0 +1,63 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2023 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. +// +// --------------------------------------------------------------------- + +// Check that we can correctly serialize a vector via boost. + +#include +#include + +#include + +#include "../tests.h" + + +int +main(int argc, char **argv) +{ + Utilities::MPI::MPI_InitFinalize mpi_context(argc, argv, 1); + MPILogInitAll mpi_log; + + const auto rank = Utilities::MPI::this_mpi_process(MPI_COMM_WORLD); + const auto n_mpi_processes = Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD); + + constexpr types::global_dof_index dofs_per_process = 5; + IndexSet owned_indices(dofs_per_process * n_mpi_processes); + owned_indices.add_range(dofs_per_process * rank, + dofs_per_process * (rank + 1)); + owned_indices.compress(); + + owned_indices.print(deallog); + + PETScWrappers::MPI::Vector petsc_vector(owned_indices, MPI_COMM_WORLD); + for (const auto dof : owned_indices) + petsc_vector(dof) = double(dof); + + petsc_vector.compress(VectorOperation::insert); + + std::stringstream stream; + boost::archive::binary_oarchive out_archive(stream); + out_archive << petsc_vector; + + PETScWrappers::MPI::Vector petsc_copy(owned_indices, MPI_COMM_WORLD); + + boost::archive::binary_iarchive in_archive(stream); + in_archive >> petsc_copy; + + petsc_copy.print(deallog.get_file_stream()); + + AssertThrow(petsc_copy == petsc_vector, ExcInternalError()); + + deallog << "OK" << std::endl; +} diff --git a/tests/petsc/serialize_vector_01.mpirun=2.output b/tests/petsc/serialize_vector_01.mpirun=2.output new file mode 100644 index 0000000000..3f56274599 --- /dev/null +++ b/tests/petsc/serialize_vector_01.mpirun=2.output @@ -0,0 +1,9 @@ + +DEAL:0::{[0,4]} +[Proc0 0-4] 0.000e+00 1.000e+00 2.000e+00 3.000e+00 4.000e+00 +DEAL:0::OK + +DEAL:1::{[5,9]} +[Proc1 5-9] 5.000e+00 6.000e+00 7.000e+00 8.000e+00 9.000e+00 +DEAL:1::OK + -- 2.39.5