/************************************************************************************* Grid physics library, www.github.com/paboyle/Grid Source file: ./lib/communicator/SharedMemory.cc Copyright (C) 2015 Author: Peter Boyle This program is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version. This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with this program; if not, write to the Free Software Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA. See the full license in the file "LICENSE" in the top level distribution directory *************************************************************************************/ /* END LEGAL */ #include #include namespace Grid { /*Construct from an MPI communicator*/ void GlobalSharedMemory::Init(Grid_MPI_Comm comm) { assert(_ShmSetup==0); WorldComm = comm; MPI_Comm_rank(WorldComm,&WorldRank); MPI_Comm_size(WorldComm,&WorldSize); // WorldComm, WorldSize, WorldRank ///////////////////////////////////////////////////////////////////// // Split into groups that can share memory ///////////////////////////////////////////////////////////////////// MPI_Comm_split_type(comm, MPI_COMM_TYPE_SHARED, 0, MPI_INFO_NULL,&WorldShmComm); MPI_Comm_rank(WorldShmComm ,&WorldShmRank); MPI_Comm_size(WorldShmComm ,&WorldShmSize); // WorldShmComm, WorldShmSize, WorldShmRank // WorldNodes WorldNodes = WorldSize/WorldShmSize; assert( (WorldNodes * WorldShmSize) == WorldSize ); // FIXME: Check all WorldShmSize are the same ? ///////////////////////////////////////////////////////////////////// // find world ranks in our SHM group (i.e. which ranks are on our node) ///////////////////////////////////////////////////////////////////// MPI_Group WorldGroup, ShmGroup; MPI_Comm_group (WorldComm, &WorldGroup); MPI_Comm_group (WorldShmComm, &ShmGroup); std::vector world_ranks(WorldSize); for(int r=0;r MyGroup; MyGroup.resize(WorldShmSize); for(int rank=0;rank()); int myleader = MyGroup[0]; std::vector leaders_1hot(WorldSize,0); std::vector leaders_group(WorldNodes,0); leaders_1hot [ myleader ] = 1; /////////////////////////////////////////////////////////////////// // global sum leaders over comm world /////////////////////////////////////////////////////////////////// int ierr=MPI_Allreduce(MPI_IN_PLACE,&leaders_1hot[0],WorldSize,MPI_INT,MPI_SUM,WorldComm); assert(ierr==0); /////////////////////////////////////////////////////////////////// // find the group leaders world rank /////////////////////////////////////////////////////////////////// int group=0; for(int l=0;l &processors,Grid_MPI_Comm & optimal_comm) { //////////////////////////////////////////////////////////////// // Assert power of two shm_size. //////////////////////////////////////////////////////////////// int log2size = -1; for(int i=0;i<=MAXLOG2RANKSPERNODE;i++){ if ( (0x1< processor_coor(ndimension); std::vector WorldDims = processors; std::vector ShmDims (ndimension,1); std::vector NodeDims (ndimension); std::vector ShmCoor (ndimension); std::vector NodeCoor (ndimension); std::vector WorldCoor(ndimension); int dim = 0; for(int l2=0;l2pw_name,WorldNode,r); shm_unlink(shm_name); int fd=shm_open(shm_name,O_RDWR|O_CREAT,0666); if ( fd < 0 ) { perror("failed shm_open"); assert(0); } ftruncate(fd, size); int mmap_flag = MAP_SHARED; #ifdef MAP_POPULATE mmap_flag |= MAP_POPULATE; #endif #ifdef MAP_HUGETLB if (flags) mmap_flag |= MAP_HUGETLB; #endif void * ptr = mmap(NULL,size, PROT_READ | PROT_WRITE, mmap_flag, fd, 0); std::cout << "Set WorldShmCommBufs["<