From 33a3fd4deaaab72bf0d7518783bd57d5c11866d4 Mon Sep 17 00:00:00 2001 From: schnellerhase <56360279+schnellerhase@users.noreply.github.com> Date: Thu, 1 Aug 2024 15:22:49 +0200 Subject: [PATCH] Introduce transfer matrix operation --- cpp/dolfinx/CMakeLists.txt | 1 + cpp/dolfinx/transfer/CMakeLists.txt | 4 + cpp/dolfinx/transfer/transfer_matrix.h | 236 +++++++++++++ cpp/test/CMakeLists.txt | 1 + cpp/test/transfer/transfer_matrix.cpp | 377 +++++++++++++++++++++ python/CMakeLists.txt | 1 + python/dolfinx/transfer.py | 23 ++ python/dolfinx/wrappers/dolfinx.cpp | 6 + python/dolfinx/wrappers/transfer.cpp | 67 ++++ python/test/unit/transfer/test_transfer.py | 22 ++ 10 files changed, 738 insertions(+) create mode 100644 cpp/dolfinx/transfer/CMakeLists.txt create mode 100644 cpp/dolfinx/transfer/transfer_matrix.h create mode 100644 cpp/test/transfer/transfer_matrix.cpp create mode 100644 python/dolfinx/transfer.py create mode 100644 python/dolfinx/wrappers/transfer.cpp create mode 100644 python/test/unit/transfer/test_transfer.py diff --git a/cpp/dolfinx/CMakeLists.txt b/cpp/dolfinx/CMakeLists.txt index 63d941fef51..d86501af4cc 100644 --- a/cpp/dolfinx/CMakeLists.txt +++ b/cpp/dolfinx/CMakeLists.txt @@ -18,6 +18,7 @@ set(DOLFINX_DIRS mesh nls refinement + transfer ) # Add source to dolfinx target, and get sets of header files diff --git a/cpp/dolfinx/transfer/CMakeLists.txt b/cpp/dolfinx/transfer/CMakeLists.txt new file mode 100644 index 00000000000..68bcddd738e --- /dev/null +++ b/cpp/dolfinx/transfer/CMakeLists.txt @@ -0,0 +1,4 @@ +set(HEADERS_transfer + ${CMAKE_CURRENT_SOURCE_DIR}/transfer_matrix.h + PARENT_SCOPE +) diff --git a/cpp/dolfinx/transfer/transfer_matrix.h b/cpp/dolfinx/transfer/transfer_matrix.h new file mode 100644 index 00000000000..add92cdc9d1 --- /dev/null +++ b/cpp/dolfinx/transfer/transfer_matrix.h @@ -0,0 +1,236 @@ +// Copyright (C) 2024 Paul T. Kühner +// +// This file is part of DOLFINX (https://www.fenicsproject.org) +// +// SPDX-License-Identifier: LGPL-3.0-or-later + +#pragma once + +#include +#include +#include +#include +#include + +#include "dolfinx/common/IndexMap.h" +#include "dolfinx/fem/FunctionSpace.h" +#include "dolfinx/la/SparsityPattern.h" +#include "dolfinx/la/utils.h" +#include "dolfinx/mesh/Mesh.h" + +namespace dolfinx::transfer +{ + +template +std::vector +inclusion_mapping(const dolfinx::mesh::Mesh& mesh_from, + const dolfinx::mesh::Mesh& mesh_to) +{ + + const common::IndexMap& im_from = *mesh_from.topology()->index_map(0); + const common::IndexMap& im_to = *mesh_to.topology()->index_map(0); + + std::vector map(im_from.size_global(), -1); + + std::span x_from = mesh_from.geometry().x(); + std::span x_to = mesh_to.geometry().x(); + + auto build_global_to_local = [&](const auto& im) + { + return [&](std::int32_t idx) + { + std::array tmp; + im.local_to_global(std::vector{idx}, tmp); + return tmp[0]; + }; + }; + + auto to_global_to = build_global_to_local(im_to); + auto to_global_from = build_global_to_local(im_from); + + for (std::int32_t i = 0; i < im_from.size_local(); i++) + { + std::ranges::subrange vertex_from(std::next(x_from.begin(), 3 * i), + std::next(x_from.begin(), 3 * (i + 1))); + for (std::int64_t j = 0; j < im_to.size_local() + im_to.num_ghosts(); j++) + { + std::ranges::subrange vertex_to(std::next(x_to.begin(), 3 * j), + std::next(x_to.begin(), 3 * (j + 1))); + + if (std::ranges::equal( + vertex_from, vertex_to, [](T a, T b) + { return std::abs(a - b) <= std::numeric_limits::epsilon(); })) + { + map[to_global_from(i)] = to_global_to(j); + break; + } + } + } + + if (dolfinx::MPI::size(mesh_to.comm()) == 1) + { + // no communication required + assert(std::ranges::all_of(map, [](auto e) { return e >= 0; })); + return map; + } + + // map holds at this point for every original local index the corresponding + // mapped global index. All other entries are still -1, but are available on + // other processes. + std::vector result(map.size(), -1); + MPI_Allreduce(map.data(), result.data(), map.size(), + dolfinx::MPI::mpi_type(), MPI_MAX, + mesh_from.comm()); + + assert(std::ranges::all_of(result, [](auto e) { return e >= 0; })); + return result; +} + +template +la::SparsityPattern +create_sparsity_pattern(const dolfinx::fem::FunctionSpace& V_from, + const dolfinx::fem::FunctionSpace& V_to, + const std::vector& inclusion_map) +{ + // TODO: P1 and which elements supported? + auto mesh_from = V_from.mesh(); + auto mesh_to = V_to.mesh(); + assert(mesh_from); + assert(mesh_to); + + MPI_Comm comm = mesh_from->comm(); + { + // Check comms equal + int result; + MPI_Comm_compare(comm, mesh_to->comm(), &result); + assert(result == MPI_CONGRUENT); + } + assert(mesh_from->topology()->dim() == mesh_to->topology()->dim()); + + auto to_v_to_f = mesh_to->topology()->connectivity(0, 1); + auto to_f_to_v = mesh_to->topology()->connectivity(1, 0); + assert(to_v_to_f); + assert(to_f_to_v); + + auto dofmap_from = V_from.dofmap(); + auto dofmap_to = V_to.dofmap(); + assert(dofmap_from); + assert(dofmap_to); + + assert(mesh_to->topology()->index_map(0)); + assert(mesh_from->topology()->index_map(0)); + const common::IndexMap& im_to = *mesh_to->topology()->index_map(0); + const common::IndexMap& im_from = *mesh_from->topology()->index_map(0); + + dolfinx::la::SparsityPattern sp( + comm, {dofmap_from->index_map, dofmap_to->index_map}, + {dofmap_from->index_map_bs(), dofmap_to->index_map_bs()}); + + assert(inclusion_map.size() == im_from.size_global()); + for (int dof_from_global = 0; dof_from_global < im_from.size_global(); + dof_from_global++) + { + std::int64_t dof_to_global = inclusion_map[dof_from_global]; + + std::vector local_dof_to_v{0}; + im_to.global_to_local(std::vector{dof_to_global}, + local_dof_to_v); + + auto local_dof_to = local_dof_to_v[0]; + + bool is_remote = (local_dof_to == -1); + bool is_ghost = local_dof_to >= im_to.size_local(); + if (is_remote || is_ghost) + continue; + + std::vector dof_from_v{0}; + im_from.global_to_local(std::vector{dof_from_global}, + dof_from_v); + + std::ranges::for_each( + to_v_to_f->links(local_dof_to), + [&](auto e) + { + sp.insert( + std::vector(to_f_to_v->links(e).size(), dof_from_v[0]), + to_f_to_v->links(e)); + }); + } + sp.finalize(); + return sp; +} + +template +void assemble_transfer_matrix(la::MatSet auto mat_set, + const dolfinx::fem::FunctionSpace& V_from, + const dolfinx::fem::FunctionSpace& V_to, + const std::vector& inclusion_map, + std::function weight) +{ + // TODO: P1 and which elements supported? + auto mesh_from = V_from.mesh(); + auto mesh_to = V_to.mesh(); + assert(mesh_from); + assert(mesh_to); + + MPI_Comm comm = mesh_from->comm(); + { + // Check comms equal + int result; + MPI_Comm_compare(comm, mesh_to->comm(), &result); + assert(result == MPI_CONGRUENT); + } + assert(mesh_from->topology()->dim() == mesh_to->topology()->dim()); + + auto to_v_to_f = mesh_to->topology()->connectivity(0, 1); + auto to_f_to_v = mesh_to->topology()->connectivity(1, 0); + assert(to_v_to_f); + assert(to_f_to_v); + + auto dofmap_from = V_from.dofmap(); + auto dofmap_to = V_to.dofmap(); + assert(dofmap_from); + assert(dofmap_to); + + assert(mesh_to->topology()->index_map(0)); + assert(mesh_from->topology()->index_map(0)); + const common::IndexMap& im_to = *mesh_to->topology()->index_map(0); + const common::IndexMap& im_from = *mesh_from->topology()->index_map(0); + + assert(inclusion_map.size() == im_from.size_global()); + + for (int dof_from_global = 0; dof_from_global < im_from.size_global(); + dof_from_global++) + { + std::int64_t dof_to_global = inclusion_map[dof_from_global]; + + std::vector local_dof_to_v{0}; + im_to.global_to_local(std::vector{dof_to_global}, + local_dof_to_v); + + auto local_dof_to = local_dof_to_v[0]; + + bool is_remote = (local_dof_to == -1); + bool is_ghost = local_dof_to >= im_to.size_local(); + if (is_remote || is_ghost) + continue; + + std::vector dof_from_v{0}; + im_from.global_to_local(std::vector{dof_from_global}, + dof_from_v); + + for (auto e : to_v_to_f->links(local_dof_to)) + { + for (auto n : to_f_to_v->links(e)) + { + // For now we only support distance <= 1 -> this should be somehow + // asserted. + std::int32_t distance = (n == local_dof_to) ? 0 : 1; + mat_set(std::vector{dof_from_v[0]}, std::vector{n}, + std::vector{weight(distance)}); + } + } + } +} + +} // namespace dolfinx::transfer \ No newline at end of file diff --git a/cpp/test/CMakeLists.txt b/cpp/test/CMakeLists.txt index 4a3dd8f3817..64315eada0e 100644 --- a/cpp/test/CMakeLists.txt +++ b/cpp/test/CMakeLists.txt @@ -47,6 +47,7 @@ add_executable( mesh/refinement/interval.cpp mesh/refinement/option.cpp mesh/refinement/rectangle.cpp + transfer/transfer_matrix.cpp ${CMAKE_CURRENT_BINARY_DIR}/poisson.c ) target_link_libraries(unittests PRIVATE Catch2::Catch2WithMain dolfinx) diff --git a/cpp/test/transfer/transfer_matrix.cpp b/cpp/test/transfer/transfer_matrix.cpp new file mode 100644 index 00000000000..3487090d01b --- /dev/null +++ b/cpp/test/transfer/transfer_matrix.cpp @@ -0,0 +1,377 @@ +// Copyright (C) 2024 Paul T. Kühner +// +// This file is part of DOLFINX (https://www.fenicsproject.org) +// +// SPDX-License-Identifier: LGPL-3.0-or-later + +#include +#include +#include +#include +#include + +#include + +#include +#include +#include +#include +#include + +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include +#include + +using namespace dolfinx; +using namespace Catch::Matchers; + +namespace +{ +template +constexpr auto EPS = std::numeric_limits::epsilon(); +} // namespace + +template +constexpr auto weight = [](std::int32_t d) -> T { return d == 0 ? 1. : .5; }; + +TEMPLATE_TEST_CASE("Transfer Matrix 1D", "[transfer_matrix]", + double) // TODO: float +{ + using T = TestType; + + auto mesh_coarse + = std::make_shared>(dolfinx::mesh::create_interval( + MPI_COMM_SELF, 2, {0.0, 1.0}, mesh::GhostMode::none)); + + auto [mesh_fine, parent_cell, parent_facet] + = refinement::refine(*mesh_coarse, std::nullopt); + + auto element = basix::create_element( + basix::element::family::P, basix::cell::type::interval, 1, + basix::element::lagrange_variant::unset, + basix::element::dpc_variant::unset, false); + + auto V_coarse = std::make_shared>( + fem::create_functionspace(mesh_coarse, element, {})); + auto V_fine + = std::make_shared>(fem::create_functionspace( + std::make_shared>(mesh_fine), element, {})); + + mesh_fine.topology()->create_connectivity(1, 0); + mesh_fine.topology()->create_connectivity(0, 1); + + std::vector inclusion_map + = transfer::inclusion_mapping(*mesh_coarse, mesh_fine); + + CHECK_THAT(inclusion_map, RangeEquals(std::vector{0, 2, 4})); + + la::SparsityPattern sp + = transfer::create_sparsity_pattern(*V_coarse, *V_fine, inclusion_map); + + la::MatrixCSR transfer_matrix(std::move(sp), la::BlockMode::compact); + transfer::assemble_transfer_matrix(transfer_matrix.mat_set_values(), + *V_coarse, *V_fine, inclusion_map, + weight); + + std::vector expected{1.0, .5, 0, 0, 0, 0, .5, 1, .5, 0, 0, 0, 0, .5, 1}; + CHECK_THAT(transfer_matrix.to_dense(), + RangeEquals(expected, [](auto a, auto b) + { return std::abs(a - b) <= EPS; })); +} + +TEMPLATE_TEST_CASE("Transfer Matrix 1D (parallel)", "[transfer_matrix]", + double) // TODO: float +{ + using T = TestType; + + int comm_size = dolfinx::MPI::size(MPI_COMM_WORLD); + int rank = dolfinx::MPI::rank(MPI_COMM_WORLD); + + // TODO: see https://github.com/FEniCS/dolfinx/issues/3358 + // auto part = mesh::create_cell_partitioner(mesh::GhostMode::shared_facet); + mesh::CellPartitionFunction part + = [&](MPI_Comm /* comm */, int /* nparts */, + const std::vector& /* cell_types */, + const std::vector>& /* cells */) + { + std::vector> data; + if (comm_size == 1) + data = {{0}, {0}, {0}, {0}}; + else if (comm_size == 2) + data = {{0}, {0}, {0}, {0}, {0, 1}, {1, 0}, {1}, {1}, {1}}; + else if (comm_size == 3) + data = {{1}, {1}, {1}, {1}, {1, 2}, {2, 1}, {2}, + {2}, {2, 0}, {0, 2}, {0}, {0}, {0}, {0}}; + else + FAIL("Test only supports <= 3 processes"); + return graph::AdjacencyList(std::move(data)); + }; + + auto mesh_coarse = std::make_shared>( + mesh::create_interval(MPI_COMM_WORLD, 5 * comm_size - 1, {0., 1.}, + mesh::GhostMode::shared_facet, part)); + + // TODO: see https://github.com/FEniCS/dolfinx/issues/3358 + // auto part_refine + // = mesh::create_cell_partitioner(mesh::GhostMode::shared_facet); + mesh::CellPartitionFunction part_refine + = [&](MPI_Comm comm, int /* nparts */, + const std::vector& /* cell_types */, + const std::vector>& /* cells */) + { + std::vector> data; + if (dolfinx::MPI::size(comm) == 1) + data = {{0}, {0}, {0}, {0}, {0}, {0}, {0}, {0}}; + else if (dolfinx::MPI::size(comm) == 2) + { + if (rank == 0) + data = {{0}, {0}, {0}, {0}, {0}, {0}, {0}, {0}, {0, 1}, {1, 0}}; + else + data = {{1}, {1}, {1}, {1}, {1}, {1}, {1}, {1}}; + } + else if (dolfinx::MPI::size(comm) == 3) + { + if (rank == 0) + data = {{2, 0}, {0, 2}, {0}, {0}, {0}, {0}, {0}, {0}, {0}, {0}}; + else if (rank == 1) + data = {{1}, {1}, {1}, {1}, {1}, {1}, {1}, {1}, {1}, {1, 2}}; + else + data = {{2, 1}, {2}, {2}, {2}, {2}, {2}, {2}, {2}}; + } + else + FAIL("Test only supports <= 3 processes"); + return graph::AdjacencyList(std::move(data)); + }; + + auto [mesh_fine, parent_cell, parent_facet] = refinement::refine( + *mesh_coarse, std::nullopt, part_refine, refinement::Option::none); + + auto element = basix::create_element( + basix::element::family::P, basix::cell::type::interval, 1, + basix::element::lagrange_variant::unset, + basix::element::dpc_variant::unset, false); + + auto V_coarse = std::make_shared>( + fem::create_functionspace(mesh_coarse, element, {})); + auto V_fine + = std::make_shared>(fem::create_functionspace( + std::make_shared>(mesh_fine), element, {})); + + mesh_fine.topology()->create_connectivity(1, 0); + mesh_fine.topology()->create_connectivity(0, 1); + + std::vector inclusion_map; + switch (comm_size) + { + case 1: + inclusion_map = {0, 2, 4, 6, 8}; + break; + case 2: + inclusion_map = {0, 2, 4, 6, 8, 10, 12, 14, 16, 18}; + break; + case 3: + inclusion_map = {0, 2, 4, 6, 8, 9, 11, 13, 15, 17, 19, 27, 25, 23, 20}; + break; + } + + CHECK_THAT(transfer::inclusion_mapping(*mesh_coarse, mesh_fine), + RangeEquals(inclusion_map)); + + la::SparsityPattern sp + = transfer::create_sparsity_pattern(*V_coarse, *V_fine, inclusion_map); + + la::MatrixCSR transfer_matrix(std::move(sp), la::BlockMode::compact); + transfer::assemble_transfer_matrix(transfer_matrix.mat_set_values(), + *V_coarse, *V_fine, inclusion_map, + weight); + + std::array, 3> expected; + if (comm_size == 1) + { + // clang-format off + expected[0] = { + /* row_0 */ 1.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_1 */ 0.0, 0.5, 1.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_2 */ 0.0, 0.0, 0.0, 0.5, 1.0, 0.5, 0.0, 0.0, 0.0, + /* row_3 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1.0, 0.5, 0.0, + /* row_4 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1.0 + }; + // clang-format on + } + else if (comm_size == 2) + { + // clang-format off + expected[0] = { + /* row_0 */ 1.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_1 */ 0.0, 0.5, 1.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_2 */ 0.0, 0.0, 0.0, 0.5, 1.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_3 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_4 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_5 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_6 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 + }; + expected[1] = { + /* row_0 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_1 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1.0, 0.5, 0.0, 0.0, 0.0, + /* row_2 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1.0, 0.5, 0.0, + /* row_3 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1.0, + /* row_4 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_5 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 + }; + // clang-format on + } + else if (comm_size == 3) + { + // clang-format off + expected[0] = { + /* row_0 */ 1.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_1 */ 0.0, 0.5, 1.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_2 */ 0.0, 0.0, 0.0, 0.5, 1.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_3 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_4 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_5 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_6 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + }; + expected[1] = { + /* row_0 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_1 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_2 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_3 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_4 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_5 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, + /* row_6 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + }; + expected[2] = { + /* row_0 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1.0, 0.5, + /* row_1 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1.0, 0.5, 0.0, 0.0, + /* row_2 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 1.0, 0.5, 0.0, 0.0, 0.0, 0.0, + /* row_3 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 1.0, 0.5, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_4 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_5 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + /* row_6 */ 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0 + }; + // clang-format on + } + else + { + CHECK(false); + } + + CHECK_THAT(transfer_matrix.to_dense(), + RangeEquals(expected[rank], [](auto a, auto b) + { return std::abs(a - b) <= EPS; })); +} + +TEMPLATE_TEST_CASE("Transfer Matrix 2D", "[transfer_matrix]", double) +{ + using T = TestType; + + auto mesh_coarse + = std::make_shared>(dolfinx::mesh::create_rectangle( + MPI_COMM_SELF, {{{0, 0}, {1, 1}}}, {1, 1}, mesh::CellType::triangle)); + mesh_coarse->topology()->create_entities(1); + + auto [mesh_fine, parent_cell, parent_facet] + = refinement::refine(*mesh_coarse, std::nullopt); + + auto element = basix::create_element( + basix::element::family::P, basix::cell::type::triangle, 1, + basix::element::lagrange_variant::unset, + basix::element::dpc_variant::unset, false); + + auto V_coarse = std::make_shared>( + fem::create_functionspace(mesh_coarse, element, {})); + auto V_fine + = std::make_shared>(fem::create_functionspace( + std::make_shared>(mesh_fine), element, {})); + + mesh_fine.topology()->create_connectivity(1, 0); + mesh_fine.topology()->create_connectivity(0, 1); + + std::vector inclusion_map + = transfer::inclusion_mapping(*mesh_coarse, mesh_fine); + CHECK_THAT(inclusion_map, RangeEquals(std::vector{4, 1, 5, 8})); + + la::SparsityPattern sp + = transfer::create_sparsity_pattern(*V_coarse, *V_fine, inclusion_map); + la::MatrixCSR transfer_matrix(std::move(sp), la::BlockMode::compact); + transfer::assemble_transfer_matrix(transfer_matrix.mat_set_values(), + *V_coarse, *V_fine, inclusion_map, + weight); + transfer_matrix.scatter_rev(); + std::vector expected{0.5, 0.0, 0.5, 0.0, 1.0, 0.0, 0.5, 0.0, 0.0, + 0.5, 1.0, 0.5, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, + 0.5, 0.0, 0.0, 0.5, 0.0, 1.0, 0.0, 0.5, 0.0, + 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.5, 1.0}; + CHECK_THAT(transfer_matrix.to_dense(), + RangeEquals(expected, [](auto a, auto b) + { return std::abs(a - b) <= EPS; })); +} + +TEMPLATE_TEST_CASE("Transfer Matrix 3D", "[transfer_matrix]", double) +{ + using T = TestType; + + auto mesh_coarse = std::make_shared>( + dolfinx::mesh::create_box(MPI_COMM_SELF, {{{0, 0, 0}, {1, 1, 1}}}, + {1, 1, 1}, mesh::CellType::tetrahedron)); + mesh_coarse->topology()->create_entities(1); + + auto [mesh_fine, parent_cell, parent_facet] + = refinement::refine(*mesh_coarse, std::nullopt); + + auto element = basix::create_element( + basix::element::family::P, basix::cell::type::tetrahedron, 1, + basix::element::lagrange_variant::unset, + basix::element::dpc_variant::unset, false); + + auto V_coarse = std::make_shared>( + fem::create_functionspace(mesh_coarse, element, {})); + auto V_fine + = std::make_shared>(fem::create_functionspace( + std::make_shared>(mesh_fine), element, {})); + + mesh_fine.topology()->create_connectivity(1, 0); + mesh_fine.topology()->create_connectivity(0, 1); + + std::vector inclusion_map + = transfer::inclusion_mapping(*mesh_coarse, mesh_fine); + CHECK_THAT(inclusion_map, RangeEquals(std::vector{ + 0, 6, 15, 25, 17, 9, 11, 22})); + + la::SparsityPattern sp + = transfer::create_sparsity_pattern(*V_coarse, *V_fine, inclusion_map); + la::MatrixCSR transfer_matrix(std::move(sp), la::BlockMode::compact); + transfer::assemble_transfer_matrix(transfer_matrix.mat_set_values(), + *V_coarse, *V_fine, inclusion_map, + weight); + + std::vector expected{ + 1.0, 0.5, 0.5, 0.5, 0.5, 0.5, 0.0, 0.5, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.5, + 0.5, 0.5, 0.0, 1.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.5, 0.0, 0.5, 0.0, 0.0, 0.0, + 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.5, 0.0, 0.0, 0.0, + 0.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, 0.5, 0.5, 1.0, 0.0, 0.0, 0.0, 0.5, 0.0, + 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.5, 0.0, 0.5, + 0.5, 1.0, 0.5, 0.0, 0.5, 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, + 0.5, 0.0, 0.5, 0.0, 0.5, 1.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, 0.5, 0.0, 0.0, + 0.0, 0.5, 0.5, 0.0, 0.0, 0.5, 0.0, 0.0, 0.5, 1.0, 0.0, 0.0, 0.0, 0.5, 0.0, + 0.0, 0.0, 0.0, 0.5, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.0, + 0.0, 0.5, 0.0, 0.0, 0.5, 0.5, 0.0, 0.0, 1.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.0, + 0.0, 0.0, 0.5, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.0, + 0.0, 0.0, 0.5, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.5, 0.5, 0.5, + 0.5, 1.0, 0.0, 0.0, 0.0, 0.5}; + CHECK_THAT(transfer_matrix.to_dense(), + RangeEquals(expected, [](auto a, auto b) + { return std::abs(a - b) <= EPS; })); +} diff --git a/python/CMakeLists.txt b/python/CMakeLists.txt index 2a9f9de56ee..90dce7571e8 100644 --- a/python/CMakeLists.txt +++ b/python/CMakeLists.txt @@ -46,6 +46,7 @@ nanobind_add_module( dolfinx/wrappers/mesh.cpp dolfinx/wrappers/petsc.cpp dolfinx/wrappers/refinement.cpp + dolfinx/wrappers/transfer.cpp ) target_compile_definitions(cpp PRIVATE cxx_std_20) diff --git a/python/dolfinx/transfer.py b/python/dolfinx/transfer.py new file mode 100644 index 00000000000..4e2a7670f3b --- /dev/null +++ b/python/dolfinx/transfer.py @@ -0,0 +1,23 @@ +import numpy as np +from numpy.typing import NDArray + +from dolfinx.cpp.la import SparsityPattern +from dolfinx.cpp.transfer import ( + assemble_transfer_matrix, +) +from dolfinx.cpp.transfer import create_sparsity_pattern as _create_sparsity_pattern +from dolfinx.cpp.transfer import inclusion_mapping as _inclusion_mapping +from dolfinx.fem import FunctionSpace +from dolfinx.mesh import Mesh + +__all__ = ["inclusion_mapping", "create_sparsity_pattern", "assemble_transfer_matrix"] + + +def inclusion_mapping(mesh_from: Mesh, mesh_to: Mesh) -> NDArray[np.int64]: + return _inclusion_mapping(mesh_from._cpp_object, mesh_to._cpp_object) + + +def create_sparsity_pattern( + V_from: FunctionSpace, V_to: FunctionSpace, inclusion_map: NDArray[np.int64] +) -> SparsityPattern: + return _create_sparsity_pattern(V_from._cpp_object, V_to._cpp_object, inclusion_map) diff --git a/python/dolfinx/wrappers/dolfinx.cpp b/python/dolfinx/wrappers/dolfinx.cpp index 9c0864f9ef1..9d2ca3e716d 100644 --- a/python/dolfinx/wrappers/dolfinx.cpp +++ b/python/dolfinx/wrappers/dolfinx.cpp @@ -25,6 +25,8 @@ void la(nb::module_& m); void mesh(nb::module_& m); void nls(nb::module_& m); void refinement(nb::module_& m); +void transfer(nb::module_& m); + } // namespace dolfinx_wrappers NB_MODULE(cpp, m) @@ -74,6 +76,10 @@ NB_MODULE(cpp, m) nb::module_ refinement = m.def_submodule("refinement", "Refinement module"); dolfinx_wrappers::refinement(refinement); + // Create transfer submodule + nb::module_ transfer = m.def_submodule("transfer", "Transfer module"); + dolfinx_wrappers::transfer(transfer); + #if defined(HAS_PETSC) && defined(HAS_PETSC4PY) // PETSc-specific wrappers nb::module_ nls = m.def_submodule("nls", "Nonlinear solver module"); diff --git a/python/dolfinx/wrappers/transfer.cpp b/python/dolfinx/wrappers/transfer.cpp new file mode 100644 index 00000000000..16ffc4cbe5a --- /dev/null +++ b/python/dolfinx/wrappers/transfer.cpp @@ -0,0 +1,67 @@ +// Copyright (C) 2024 Paul T. Kühner +// +// This file is part of DOLFINX (https://www.fenicsproject.org) +// +// SPDX-License-Identifier: LGPL-3.0-or-later + +#include +#include + +#include +#include + +#include +#include + +#include "array.h" + +namespace nb = nanobind; + +namespace dolfinx_wrappers +{ + +void transfer(nb::module_& m) +{ + m.def( + "inclusion_mapping", + [](const dolfinx::mesh::Mesh& mesh_from, + const dolfinx::mesh::Mesh& mesh_to) + { + std::vector map + = dolfinx::transfer::inclusion_mapping(mesh_from, mesh_to); + return dolfinx_wrappers::as_nbarray(map); + }, + nb::arg("mesh_from"), nb::arg("mesh_to"), + "Computes inclusion mapping between two meshes"); + + m.def( + "create_sparsity_pattern", + [](const dolfinx::fem::FunctionSpace& V_from, + const dolfinx::fem::FunctionSpace& V_to, + nb::ndarray& inclusion_map) + { + // TODO: cahnge to accepting range; + auto vec = std::vector(inclusion_map.data(), + inclusion_map.data() + inclusion_map.size()); + return dolfinx::transfer::create_sparsity_pattern(V_from, V_to, + vec); + }, + nb::arg("V_from"), nb::arg("V_to"), nb::arg("inclusion_map")); + + m.def( + "assemble_transfer_matrix", + [](dolfinx::la::MatrixCSR& A, + const dolfinx::fem::FunctionSpace& V_from, + const dolfinx::fem::FunctionSpace& V_to, + const std::vector& inclusion_map, + std::function weight) + { + dolfinx::transfer::assemble_transfer_matrix( + A.mat_set_values(), V_from, V_to, inclusion_map, weight); + }, + nb::arg("A"), nb::arg("V_from"), nb::arg("V_to"), + nb::arg("inclusion_map"), nb::arg("weight"), + "Assembles a transfer matrix between two function spaces."); +} + +} // namespace dolfinx_wrappers diff --git a/python/test/unit/transfer/test_transfer.py b/python/test/unit/transfer/test_transfer.py new file mode 100644 index 00000000000..93c260442da --- /dev/null +++ b/python/test/unit/transfer/test_transfer.py @@ -0,0 +1,22 @@ +from mpi4py import MPI + +import pytest + +from dolfinx.fem import functionspace +from dolfinx.mesh import GhostMode, create_interval, refine +from dolfinx.transfer import create_sparsity_pattern, inclusion_mapping + + +@pytest.mark.parametrize( + "ghost_mode", [GhostMode.none, GhostMode.shared_vertex, GhostMode.shared_facet] +) +def test_1d(ghost_mode): + mesh = create_interval(MPI.COMM_WORLD, 10, (0, 1), ghost_mode=ghost_mode) + mesh_fine, _, _ = refine(mesh) + inclusion_map = inclusion_mapping(mesh, mesh_fine) + V = functionspace(mesh, ("Lagrange", 1, (1,))) + V_fine = functionspace(mesh_fine, ("Lagrange", 1, (1,))) + V_fine.mesh.topology.create_connectivity(1, 0) + V_fine.mesh.topology.create_connectivity(0, 1) + create_sparsity_pattern(V, V_fine, inclusion_map) + # continue with assembly of matrix