mirror of
https://github.com/cosmo-sims/monofonIC.git
synced 2024-09-19 17:03:45 +02:00
fixed mpi_init_thread to allow for MPI_THREAD_MULTIPLE allowing multithread receives
This commit is contained in:
parent
0420726270
commit
d8cd286f87
2 changed files with 21 additions and 15 deletions
|
@ -340,7 +340,7 @@ void Grid_FFT<data_t, bdistributed>::FourierInterpolateCopyTo( grid_fft_t &grid_
|
|||
MPI_Allgather(&grid_to.local_1_start_, 1, MPI_LONG_LONG, &offsets_recv[0], 1,
|
||||
MPI_LONG_LONG, MPI_COMM_WORLD);
|
||||
|
||||
MPI_Datatype datatype =
|
||||
const MPI_Datatype datatype =
|
||||
(typeid(data_t) == typeid(float)) ? MPI_C_FLOAT_COMPLEX
|
||||
: (typeid(data_t) == typeid(double)) ? MPI_C_DOUBLE_COMPLEX
|
||||
: (typeid(data_t) == typeid(long double)) ? MPI_C_LONG_DOUBLE_COMPLEX
|
||||
|
@ -349,11 +349,11 @@ void Grid_FFT<data_t, bdistributed>::FourierInterpolateCopyTo( grid_fft_t &grid_
|
|||
const size_t slicesz = grid_from.size(1) * grid_from.size(3);
|
||||
|
||||
// determine effective Nyquist modes representable by both fields and their locations in array
|
||||
size_t fny0_left = std::min(grid_from.n_[1] / 2, grid_to.n_[1] / 2);
|
||||
size_t fny0_right = std::max(grid_from.n_[1] - grid_to.n_[1] / 2, grid_from.n_[1] / 2);
|
||||
size_t fny1_left = std::min(grid_from.n_[0] / 2, grid_to.n_[0] / 2);
|
||||
size_t fny1_right = std::max(grid_from.n_[0] - grid_to.n_[0] / 2, grid_from.n_[0] / 2);
|
||||
size_t fny2_left = std::min(grid_from.n_[2] / 2, grid_to.n_[2] / 2);
|
||||
const size_t fny0_left = std::min(grid_from.n_[1] / 2, grid_to.n_[1] / 2);
|
||||
const size_t fny0_right = std::max(grid_from.n_[1] - grid_to.n_[1] / 2, grid_from.n_[1] / 2);
|
||||
const size_t fny1_left = std::min(grid_from.n_[0] / 2, grid_to.n_[0] / 2);
|
||||
const size_t fny1_right = std::max(grid_from.n_[0] - grid_to.n_[0] / 2, grid_from.n_[0] / 2);
|
||||
const size_t fny2_left = std::min(grid_from.n_[2] / 2, grid_to.n_[2] / 2);
|
||||
// size_t fny2_right = std::max(grid_from.n_[2] - grid_to.n_[2] / 2, grid_from.n_[2] / 2);
|
||||
|
||||
const size_t fny0_left_recv = fny0_left;
|
||||
|
@ -392,15 +392,17 @@ void Grid_FFT<data_t, bdistributed>::FourierInterpolateCopyTo( grid_fft_t &grid_
|
|||
}
|
||||
}
|
||||
|
||||
MPI_Barrier( MPI_COMM_WORLD );
|
||||
|
||||
//--- receive data ------------------------------------------------------------
|
||||
#pragma omp parallel
|
||||
{
|
||||
MPI_Status status;
|
||||
ccomplex_t * recvbuf = new ccomplex_t[ slicesz ];
|
||||
|
||||
#pragma omp for
|
||||
#pragma omp for schedule(dynamic)
|
||||
for( size_t i=0; i<grid_to.size(0); ++i )
|
||||
{
|
||||
{
|
||||
size_t iglobal_recv = i + offsets_recv[CONFIG::MPI_task_rank];
|
||||
|
||||
if (iglobal_recv < fny0_left_recv || iglobal_recv > fny0_right_recv)
|
||||
|
@ -409,10 +411,13 @@ void Grid_FFT<data_t, bdistributed>::FourierInterpolateCopyTo( grid_fft_t &grid_
|
|||
|
||||
int recvfrom = get_task(iglobal_send, offsets_send, CONFIG::MPI_task_size);
|
||||
|
||||
// receive data slice and check for MPI errors when in debug mode
|
||||
status.MPI_ERROR = MPI_SUCCESS;
|
||||
MPI_Recv(&recvbuf[0], (int)slicesz, datatype, recvfrom, (int)iglobal_recv, MPI_COMM_WORLD, &status);
|
||||
assert(status.MPI_ERROR == MPI_SUCCESS);
|
||||
//#pragma omp critical // need critical region here if we do "MPI_THREAD_FUNNELED",
|
||||
{
|
||||
// receive data slice and check for MPI errors when in debug mode
|
||||
status.MPI_ERROR = MPI_SUCCESS;
|
||||
MPI_Recv(&recvbuf[0], (int)slicesz, datatype, recvfrom, (int)iglobal_recv, MPI_COMM_WORLD, &status);
|
||||
assert(status.MPI_ERROR == MPI_SUCCESS);
|
||||
}
|
||||
|
||||
// copy data slice into new grid, avoiding modes that do not exist in both
|
||||
for( size_t j=0; j<grid_to.size(1); ++j )
|
||||
|
|
|
@ -74,8 +74,9 @@ int main( int argc, char** argv )
|
|||
//------------------------------------------------------------------------------
|
||||
|
||||
#if defined(USE_MPI)
|
||||
MPI_Init_thread(&argc, &argv, MPI_THREAD_FUNNELED, &CONFIG::MPI_thread_support);
|
||||
CONFIG::MPI_threads_ok = CONFIG::MPI_thread_support >= MPI_THREAD_FUNNELED;
|
||||
int thread_wanted = MPI_THREAD_MULTIPLE; // MPI_THREAD_FUNNELED
|
||||
MPI_Init_thread(&argc, &argv, thread_wanted, &CONFIG::MPI_thread_support);
|
||||
CONFIG::MPI_threads_ok = CONFIG::MPI_thread_support >= thread_wanted;
|
||||
MPI_Comm_rank(MPI_COMM_WORLD, &CONFIG::MPI_task_rank);
|
||||
MPI_Comm_size(MPI_COMM_WORLD, &CONFIG::MPI_task_size);
|
||||
CONFIG::MPI_ok = true;
|
||||
|
@ -168,7 +169,7 @@ int main( int argc, char** argv )
|
|||
omp_set_num_threads(CONFIG::num_threads);
|
||||
#endif
|
||||
|
||||
// std::feclearexcept(FE_ALL_EXCEPT);
|
||||
std::feclearexcept(FE_ALL_EXCEPT);
|
||||
|
||||
//------------------------------------------------------------------------------
|
||||
// Write code configuration to screen
|
||||
|
|
Loading…
Reference in a new issue