diff --git a/Tests/MathLib/TestGlobalVectorInterface.cpp b/Tests/MathLib/TestGlobalVectorInterface.cpp index 755cf570945ab224380a7353f6ba31bd6633afac..dca634e5a221c3aa4ac13fb98e7faad4e4912ab9 100644 --- a/Tests/MathLib/TestGlobalVectorInterface.cpp +++ b/Tests/MathLib/TestGlobalVectorInterface.cpp @@ -15,6 +15,10 @@ #include <gtest/gtest.h> +#include <cmath> +#include <numeric> +#include <ranges> + #include "MathLib/LinAlg/LinAlg.h" #include "Tests/TestTools.h" @@ -78,115 +82,100 @@ void checkGlobalVectorInterface() } #ifdef USE_PETSC -template <class T_VECTOR> -void checkGlobalVectorInterfacePETSc() +double getEuclideanNorm(std::vector<double> const& v) { - int msize; - MPI_Comm_size(PETSC_COMM_WORLD, &msize); + // Create a view that squares each element + auto squared_view = + v | std::views::transform([](double x) { return x * x; }); - ASSERT_EQ(3u, msize); - - // ----------------------------------------------------------------- - // PETSc determined partitioning - T_VECTOR x(16); + return std::sqrt( + std::accumulate(squared_view.begin(), squared_view.end(), 0.0)); +} - ASSERT_EQ(16u, x.size()); +template <class T_VECTOR> +void checkPETScVectorNoExplictGhostID(T_VECTOR& x, int const msize) +{ ASSERT_EQ(x.getRangeEnd() - x.getRangeBegin(), x.getLocalSize()); - const int r0 = x.getRangeBegin(); - // x.get(0) is expensive, only get local value. Use it for test purpose + // Note id0_rank may vary because of the PETSc determined partitioning + const int id0_rank = x.getRangeBegin(); x.setLocalAccessibleVector(); - ASSERT_EQ(.0, x.get(r0)); + ASSERT_EQ(.0, x.get(id0_rank)); - set(x, 10.); + set(x, 10.); // x:=10.0 + x.finalizeAssembly(); // Value of x is not copied to y const bool deep_copy = false; - T_VECTOR y(x, deep_copy); + T_VECTOR y(x, deep_copy); // y:= 0.0 y.setLocalAccessibleVector(); - ASSERT_EQ(0, y.get(r0)); + ASSERT_EQ(0, y.get(id0_rank)); - set(y, 10.0); + set(y, 10.0); // y:=10.0 y.setLocalAccessibleVector(); - ASSERT_EQ(10, y.get(r0)); + ASSERT_EQ(10, y.get(id0_rank)); // y += x - axpy(y, 1., x); + axpy(y, 1., x); // y = x + y := 20 y.setLocalAccessibleVector(); - ASSERT_EQ(20, y.get(r0)); - ASSERT_EQ(80., norm2(y)); + ASSERT_EQ(20, y.get(id0_rank)); // y -= x - axpy(y, -1., x); + axpy(y, -1., x); // y = y - x : = 10.0 y.setLocalAccessibleVector(); - ASSERT_EQ(10, y.get(r0)); - ASSERT_EQ(40., norm2(y)); + ASSERT_EQ(10, y.get(id0_rank)); + // Set values by each rank std::vector<double> local_vec(2, 10.0); std::vector<GlobalIndexType> vec_pos(2); - vec_pos[0] = r0; // any index in [0,15] - vec_pos[1] = r0 + 1; // any index in [0,15] - + vec_pos[0] = id0_rank; + int const id0_pls_1 = id0_rank + 1; + // If id0_pls_1 is out of the current rank, set its negative value to + // vec_pos[1] to avoid adding value. + vec_pos[1] = id0_pls_1 < y.getRangeEnd() ? id0_pls_1 : -id0_pls_1; y.add(vec_pos, local_vec); + y.finalizeAssembly(); - double normy = std::sqrt(6.0 * 400 + 10.0 * 100); + int const x_size = x.size(); + // Count the number of the change elements of y by y.add + std::vector<int> all_vec_pos(msize * 2); + MPI_Allgather(vec_pos.data(), 2, MPI_INT, all_vec_pos.data(), 2, MPI_INT, + PETSC_COMM_WORLD); - ASSERT_NEAR(0.0, normy - norm2(y), 1.e-10); + std::vector<double> x_raw(x_size); + x.getGlobalVector(x_raw); + std::ranges::for_each( + all_vec_pos, + [&x_raw](int id) + { + if (id >= 0 && id < static_cast<int>(x_raw.size())) + { + x_raw[id] = 20.0; + } + }); - std::vector<double> x0(16); - double z[] = { - 2.0000000000000000e+01, 2.0000000000000000e+01, 1.0000000000000000e+01, - 1.0000000000000000e+01, 1.0000000000000000e+01, 1.0000000000000000e+01, - 2.0000000000000000e+01, 2.0000000000000000e+01, 1.0000000000000000e+01, - 1.0000000000000000e+01, 1.0000000000000000e+01, 2.0000000000000000e+01, - 2.0000000000000000e+01, 1.0000000000000000e+01, 1.0000000000000000e+01, - 1.0000000000000000e+01}; + std::vector<double> y_raw(x_size); - y.getGlobalVector(x0); + y.getGlobalVector(y_raw); - ASSERT_ARRAY_NEAR(z, x0, 16, 1e-10); - - // ----------------------------------------------------------------- - // User determined partitioning - const bool is_global_size = false; - T_VECTOR x_fixed_p(2, is_global_size); + ASSERT_ARRAY_NEAR(x_raw, y_raw, x_size, 1e-10); - ASSERT_EQ(6u, x_fixed_p.size()); + ASSERT_NEAR(getEuclideanNorm(y_raw), norm2(y), 1e-10); - int mrank; - MPI_Comm_rank(PETSC_COMM_WORLD, &mrank); - - ASSERT_EQ(2 * mrank, x_fixed_p.getRangeBegin()); - ASSERT_EQ(2 * mrank + 2, x_fixed_p.getRangeEnd()); - - vec_pos[0] = 2 * mrank; - vec_pos[1] = vec_pos[0] + 1; - local_vec[0] = 1.; - local_vec[1] = 2.; - for (unsigned i = 0; i < 3; i++) - { - const unsigned j = 2 * i; - z[j] = 1.0; - z[j + 1] = 2.0; - } - x_fixed_p.set(vec_pos, local_vec); - x_fixed_p.getGlobalVector(x0); - - ASSERT_ARRAY_NEAR(z, x0, 6, 1e-10); - - // check local array - std::vector<double> loc_v; - x_fixed_p.copyValues(loc_v); - z[0] = 1.0; - z[1] = 2.0; + // Deep copy + T_VECTOR x_deep_copied(y); + x_deep_copied.getGlobalVector(y_raw); + ASSERT_ARRAY_NEAR(x_raw, y_raw, x_size, 1e-10); +} - ASSERT_ARRAY_NEAR(z, loc_v, 2, 1e-10); +template <class T_VECTOR> +void checkPETScVectorExplictGhostID() +{ + int msize; + MPI_Comm_size(PETSC_COMM_WORLD, &msize); - // Deep copy - MathLib::finalizeVectorAssembly(x_fixed_p); - T_VECTOR x_deep_copied(x_fixed_p); - ASSERT_NEAR(sqrt(3.0 * 5), norm2(x_deep_copied), 1.e-10); + ASSERT_EQ(3u, msize); // ----------------------------------------------------------------- // Vector with ghost entries @@ -210,16 +199,24 @@ void checkGlobalVectorInterfacePETSc() The above ghost entry embedded vector is realized by the following test. */ + int mrank; + MPI_Comm_rank(PETSC_COMM_WORLD, &mrank); std::size_t local_vec_size = 4; if (mrank == 1) + { local_vec_size = 5; + } else if (mrank == 2) + { local_vec_size = 3; + } std::vector<GlobalIndexType> non_ghost_ids(local_vec_size); std::vector<double> non_ghost_vals(local_vec_size); std::size_t nghosts = 3; - if (mrank) + if (mrank > 0) + { nghosts = 2; + } std::vector<GlobalIndexType> ghost_ids(nghosts); std::vector<double> expected; switch (mrank) @@ -243,7 +240,7 @@ void checkGlobalVectorInterfacePETSc() expected = {9., 10., 11., 3., 5.}; break; } - T_VECTOR x_with_ghosts(local_vec_size, ghost_ids, is_global_size); + T_VECTOR x_with_ghosts(local_vec_size, ghost_ids, false /*is_global_size*/); x_with_ghosts.set(non_ghost_ids, non_ghost_vals); MathLib::finalizeVectorAssembly(x_with_ghosts); @@ -263,9 +260,40 @@ void checkGlobalVectorInterfacePETSc() //-------------------------------------------- #if defined(USE_PETSC) -TEST(MPITest_Math, CheckInterface_PETScVector) +TEST(MPITest_Math, PETScVectorPatitionedAutomatically) +{ + int msize; + MPI_Comm_size(PETSC_COMM_WORLD, &msize); + + MathLib::PETScVector x(16); // x:=0.0 + ASSERT_LE(msize, 16u); + + ASSERT_EQ(16u, x.size()); + checkPETScVectorNoExplictGhostID<MathLib::PETScVector>(x, msize); +} + +TEST(MPITest_Math, PETScVectorFixedPartition) +{ + int msize; + MPI_Comm_size(PETSC_COMM_WORLD, &msize); + + // Each partition has two elements. + int const num_elements_per_rank = 2; + MathLib::PETScVector x_fixed_p(num_elements_per_rank, + false /*is_global_size*/); + + ASSERT_EQ(num_elements_per_rank * msize, x_fixed_p.size()); + + int mrank; + MPI_Comm_rank(PETSC_COMM_WORLD, &mrank); + ASSERT_EQ(2 * mrank, x_fixed_p.getRangeBegin()); + ASSERT_EQ(2 * mrank + 2, x_fixed_p.getRangeEnd()); + checkPETScVectorNoExplictGhostID<MathLib::PETScVector>(x_fixed_p, msize); +} + +TEST(MPITest_Math, CheckPETScVectorExplictGhostID) { - checkGlobalVectorInterfacePETSc<MathLib::PETScVector>(); + checkPETScVectorExplictGhostID<MathLib::PETScVector>(); } #else TEST(Math, CheckInterface_EigenVector)