source: XIOS/dev/branch_openmp/extern/ep_dev/ep_bcast.cpp @ 1527

Last change on this file since 1527 was 1527, checked in by yushan, 6 years ago

save dev

File size: 2.2 KB
RevLine 
[1381]1/*!
2   \file ep_bcast.cpp
3   \since 2 may 2016
4
5   \brief Definitions of MPI collective function: MPI_Bcast
6 */
7
8#include "ep_lib.hpp"
9#include <mpi.h>
10#include "ep_declaration.hpp"
11#include "ep_mpi.hpp"
12
13using namespace std;
14
15namespace ep_lib
16{
17
18  int MPI_Bcast_local(void *buffer, int count, MPI_Datatype datatype, int local_root, MPI_Comm comm)
19  {
20    assert(valid_type(datatype));
21
[1500]22    int ep_rank_loc = comm->ep_comm_ptr->size_rank_info[1].first;
[1381]23
24    ::MPI_Aint datasize, lb;
25    ::MPI_Type_get_extent(to_mpi_type(datatype), &lb, &datasize);
26   
27
28    if(ep_rank_loc == local_root)
29    {
[1500]30      comm->my_buffer->void_buffer[local_root] = buffer;
[1381]31    }
32
33    MPI_Barrier_local(comm);
34
35    if(ep_rank_loc != local_root)
36    {
37      #pragma omp critical (_bcast)     
[1500]38      memcpy(buffer, comm->my_buffer->void_buffer[local_root], datasize * count);
[1381]39    }
40
41    MPI_Barrier_local(comm);
42  }
43
44  int MPI_Bcast(void *buffer, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
45  {
46
[1527]47    if(!comm->is_ep) return ::MPI_Bcast(buffer, count, to_mpi_type(datatype), root, to_mpi_comm(comm->mpi_comm));
48    if(comm->is_intercomm) return MPI_Bcast_intercomm(buffer, count, datatype, root, comm);
[1381]49
50
[1500]51    int ep_rank = comm->ep_comm_ptr->size_rank_info[0].first;
52    int ep_rank_loc = comm->ep_comm_ptr->size_rank_info[1].first;
53    int mpi_rank = comm->ep_comm_ptr->size_rank_info[2].first;
[1381]54
[1503]55    int root_mpi_rank = comm->ep_rank_map->at(root).second;
56    int root_ep_rank_loc = comm->ep_rank_map->at(root).first;
[1381]57
[1527]58    //printf("ep_rank = %d, root_mpi_rank = %d, root_ep_rank_loc = %d\n", ep_rank, root_mpi_rank, root_ep_rank_loc);
[1381]59
[1527]60
[1381]61    if((ep_rank_loc==0 && mpi_rank != root_mpi_rank ) || ep_rank == root)
62    {
[1500]63      ::MPI_Bcast(buffer, count, to_mpi_type(datatype), root_mpi_rank, to_mpi_comm(comm->mpi_comm));
[1381]64    }
65
66    if(mpi_rank == root_mpi_rank) MPI_Bcast_local(buffer, count, datatype, root_ep_rank_loc, comm);
67    else                          MPI_Bcast_local(buffer, count, datatype, 0, comm);
68
69  }
70
71
[1527]72  int MPI_Bcast_intercomm(void *buffer, int count, MPI_Datatype datatype, int root, MPI_Comm comm)
73  {
74    printf("MPI_Bcast_intercomm not yet implemented\n");
75    MPI_Abort(comm, 0);
76  }
[1381]77
78
79}
Note: See TracBrowser for help on using the repository browser.