OpenPFC  0.1.4
Phase Field Crystal simulation framework
Loading...
Searching...
No Matches
exchange.hpp
Go to the documentation of this file.
1// SPDX-FileCopyrightText: 2025 VTT Technical Research Centre of Finland Ltd
2// SPDX-License-Identifier: AGPL-3.0-or-later
3
44#pragma once
45
46#include <mpi.h>
47#include <vector>
48
51
52#if defined(OpenPFC_ENABLE_CUDA)
53#include <cuda_runtime.h>
54#endif
55
56namespace pfc {
57namespace exchange {
58namespace detail {
59
63template <typename T> MPI_Datatype get_mpi_type();
64
65template <> inline MPI_Datatype get_mpi_type<double>() { return MPI_DOUBLE; }
66
67template <> inline MPI_Datatype get_mpi_type<float>() { return MPI_FLOAT; }
68
69template <> inline MPI_Datatype get_mpi_type<int>() { return MPI_INT; }
70
71template <> inline MPI_Datatype get_mpi_type<size_t>() {
72 return MPI_UNSIGNED_LONG_LONG;
73}
74
75} // namespace detail
76
89template <typename BackendTag, typename T>
91 int receiver_rank, MPI_Comm comm, int tag = 0) {
92 int my_rank;
93 MPI_Comm_rank(comm, &my_rank);
94
95 if (my_rank != sender_rank) {
96 return; // Not the sender
97 }
98
99 size_t size = sparse_vector.size();
100
101 // Send size first
103
104 if (size == 0) {
105 return;
106 }
107
108 // Get indices and data (may need to copy from device for CUDA)
109 std::vector<size_t> indices;
110 std::vector<T> data;
111
112 if constexpr (std::is_same_v<BackendTag, backend::CpuTag>) {
113 // CPU: Direct access
114 indices.resize(size);
115 data.resize(size);
116 std::copy(sparse_vector.indices().data(), sparse_vector.indices().data() + size,
117 indices.begin());
118 std::copy(sparse_vector.data().data(), sparse_vector.data().data() + size,
119 data.begin());
120 }
121#if defined(OpenPFC_ENABLE_CUDA)
122 else if constexpr (std::is_same_v<BackendTag, backend::CudaTag>) {
123 // CUDA: Copy to host first
124 indices.resize(size);
125 data.resize(size);
126 cudaMemcpy(indices.data(), sparse_vector.indices().data(), size * sizeof(size_t),
128 cudaMemcpy(data.data(), sparse_vector.data().data(), size * sizeof(T),
130 }
131#endif
132
133 // Send indices
134 MPI_Send(indices.data(), static_cast<int>(size), MPI_UNSIGNED_LONG_LONG,
135 receiver_rank, tag + 1, comm);
136
137 // Send data
138 MPI_Datatype mpi_type = detail::get_mpi_type<T>();
139 MPI_Send(data.data(), static_cast<int>(size), mpi_type, receiver_rank, tag + 2,
140 comm);
141}
142
154template <typename BackendTag, typename T>
156 int receiver_rank, MPI_Comm comm, int tag = 0) {
157 int my_rank;
158 MPI_Comm_rank(comm, &my_rank);
159
160 if (my_rank != receiver_rank) {
161 return; // Not the receiver
162 }
163
164 // Receive size first
165 size_t size;
168
169 if (size == 0) {
171 return;
172 }
173
174 // Receive indices and data
175 std::vector<size_t> indices(size);
176 std::vector<T> data(size);
177
178 int count = static_cast<int>(size);
179 MPI_Recv(indices.data(), count, MPI_UNSIGNED_LONG_LONG, sender_rank, tag + 1, comm,
181
182 MPI_Datatype mpi_type = detail::get_mpi_type<T>();
183 MPI_Recv(data.data(), count, mpi_type, sender_rank, tag + 2, comm,
185
186 // Create SparseVector from received data (will sort indices)
188}
189
202template <typename BackendTag, typename T>
204 int sender_rank, int receiver_rank, MPI_Comm comm, int tag = 0) {
205 int my_rank;
206 MPI_Comm_rank(comm, &my_rank);
207
208 if (my_rank != sender_rank) {
209 return; // Not the sender
210 }
211
212 size_t size = sparse_vector.size();
213
214 if (size == 0) {
215 return;
216 }
217
218 // Get data (may need to copy from device for CUDA)
219 std::vector<T> data;
220
221 if constexpr (std::is_same_v<BackendTag, backend::CpuTag>) {
222 // CPU: Direct access
223 data.resize(size);
224 std::copy(sparse_vector.data().data(), sparse_vector.data().data() + size,
225 data.begin());
226 }
227#if defined(OpenPFC_ENABLE_CUDA)
228 else if constexpr (std::is_same_v<BackendTag, backend::CudaTag>) {
229 // CUDA: Copy to host first
230 data.resize(size);
231 cudaMemcpy(data.data(), sparse_vector.data().data(), size * sizeof(T),
233 }
234#endif
235
236 // Send data only
237 MPI_Datatype mpi_type = detail::get_mpi_type<T>();
238 MPI_Send(data.data(), static_cast<int>(size), mpi_type, receiver_rank, tag, comm);
239}
240
253template <typename BackendTag, typename T>
255 int receiver_rank, MPI_Comm comm, int tag = 0) {
256 int my_rank;
257 MPI_Comm_rank(comm, &my_rank);
258
259 if (my_rank != receiver_rank) {
260 return; // Not the receiver
261 }
262
263 size_t size = sparse_vector.size();
264
265 if (size == 0) {
266 return;
267 }
268
269 // Receive data
270 std::vector<T> data(size);
271
272 MPI_Datatype mpi_type = exchange::detail::get_mpi_type<T>();
273 int count = static_cast<int>(size);
274 MPI_Recv(data.data(), count, mpi_type, sender_rank, tag, comm, MPI_STATUS_IGNORE);
275
276 // Copy to device if needed
277 if constexpr (std::is_same_v<BackendTag, backend::CpuTag>) {
278 // CPU: Direct copy
279 std::copy(data.begin(), data.end(), sparse_vector.data().data());
280 }
281#if defined(OpenPFC_ENABLE_CUDA)
282 else if constexpr (std::is_same_v<BackendTag, backend::CudaTag>) {
283 // CUDA: Copy to device
284 cudaMemcpy(sparse_vector.data().data(), data.data(), size * sizeof(T),
286 }
287#endif
288}
289
290} // namespace exchange
291} // namespace pfc
Backend tags for compile-time backend selection.
Sparse vector for indexed data views.
Definition sparse_vector.hpp:66
void receive(core::SparseVector< BackendTag, T > &sparse_vector, int sender_rank, int receiver_rank, MPI_Comm comm, int tag=0)
Receive SparseVector (both indices and data) - setup phase.
Definition exchange.hpp:155
void send_data(const core::SparseVector< BackendTag, T > &sparse_vector, int sender_rank, int receiver_rank, MPI_Comm comm, int tag=0)
Send only data values (indices already known) - runtime phase.
Definition exchange.hpp:203
void send(const core::SparseVector< BackendTag, T > &sparse_vector, int sender_rank, int receiver_rank, MPI_Comm comm, int tag=0)
Send SparseVector (both indices and data) - setup phase.
Definition exchange.hpp:90
MPI_Datatype get_mpi_type()
Get MPI datatype for type T.
void receive_data(core::SparseVector< BackendTag, T > &sparse_vector, int sender_rank, int receiver_rank, MPI_Comm comm, int tag=0)
Receive only data values (indices already known) - runtime phase.
Definition exchange.hpp:254
Sparse vector for halo exchange and indexed data views.
Represents the global simulation domain (the "world").
Definition world.hpp:91