OpenPFC  0.1.4
Phase Field Crystal simulation framework
Loading...
Searching...
No Matches
exchange.hpp File Reference

MPI exchange operations for SparseVector. More...

#include <mpi.h>
#include <vector>
#include <openpfc/core/backend_tags.hpp>
#include <openpfc/core/sparse_vector.hpp>
Include dependency graph for exchange.hpp:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

template<typename T >
MPI_Datatype pfc::exchange::detail::get_mpi_type ()
 Get MPI datatype for type T.
 
template<>
MPI_Datatype pfc::exchange::detail::get_mpi_type< double > ()
 
template<>
MPI_Datatype pfc::exchange::detail::get_mpi_type< float > ()
 
template<>
MPI_Datatype pfc::exchange::detail::get_mpi_type< int > ()
 
template<>
MPI_Datatype pfc::exchange::detail::get_mpi_type< size_t > ()
 
template<typename BackendTag , typename T >
void pfc::exchange::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.
 
template<typename BackendTag , typename T >
void pfc::exchange::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.
 
template<typename BackendTag , typename T >
void pfc::exchange::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.
 
template<typename BackendTag , typename T >
void pfc::exchange::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.
 

Detailed Description

MPI exchange operations for SparseVector.

Provides two-phase exchange for SparseVector:

  1. Setup phase (expensive, done once):
    • send(sparse_vec, sender_rank, receiver_rank, comm): Sends both indices and data
    • receive(sparse_vec, sender_rank, receiver_rank, comm): Receives both indices and data
  2. Runtime phase (cheap, done every step):
    • send_data(sparse_vec, sender_rank, receiver_rank, comm): Sends only data values
    • receive_data(sparse_vec, sender_rank, receiver_rank, comm): Receives only data values

Key optimization: Indices are exchanged once, then only values are transferred.

// Setup: Exchange indices once
auto sparse = sparsevector::create<double>({0, 2, 4});
exchange::send(sparse, my_rank, neighbor_rank, MPI_COMM_WORLD);
exchange::receive(sparse, neighbor_rank, my_rank, MPI_COMM_WORLD);
// Runtime: Exchange only values (repeatedly)
for (int step = 0; step < 1000; ++step) {
exchange::send_data(sparse, my_rank, neighbor_rank, MPI_COMM_WORLD);
exchange::receive_data(sparse, neighbor_rank, my_rank, MPI_COMM_WORLD);
}
See also
core/sparse_vector.hpp for SparseVector definition
Author
OpenPFC Development Team
Date
2025

Function Documentation

◆ receive()

template<typename BackendTag , typename T >
void pfc::exchange::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.

Receives both indices and data from sender. Used during initialization.

Parameters
sparse_vectorSparseVector to receive into (will be resized)
sender_rankMPI rank of sender
receiver_rankMPI rank of receiver
commMPI communicator
tagMPI message tag (default: 0)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ receive_data()

template<typename BackendTag , typename T >
void pfc::exchange::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.

Receives only the data values, assuming indices are already set. Much cheaper than full receive.

Parameters
sparse_vectorSparseVector to receive data into
sender_rankMPI rank of sender
receiver_rankMPI rank of receiver
commMPI communicator
tagMPI message tag (default: 0)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ send()

template<typename BackendTag , typename T >
void pfc::exchange::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.

Sends both indices and data to receiver. Used during initialization to establish communication pattern.

Parameters
sparse_vectorSparseVector to send
sender_rankMPI rank of sender
receiver_rankMPI rank of receiver
commMPI communicator
tagMPI message tag (default: 0)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ send_data()

template<typename BackendTag , typename T >
void pfc::exchange::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.

Sends only the data values, assuming receiver already knows the indices. Much cheaper than full send.

Parameters
sparse_vectorSparseVector to send data from
sender_rankMPI rank of sender
receiver_rankMPI rank of receiver
commMPI communicator
tagMPI message tag (default: 0)
Here is the call graph for this function:
Here is the caller graph for this function: