Grid 0.7.0
SharedMemory.cc
Go to the documentation of this file.
1/*************************************************************************************
2
3 Grid physics library, www.github.com/paboyle/Grid
4
5 Source file: ./lib/communicator/SharedMemory.cc
6
7 Copyright (C) 2015
8
9Author: Peter Boyle <paboyle@ph.ed.ac.uk>
10
11 This program is free software; you can redistribute it and/or modify
12 it under the terms of the GNU General Public License as published by
13 the Free Software Foundation; either version 2 of the License, or
14 (at your option) any later version.
15
16 This program is distributed in the hope that it will be useful,
17 but WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19 GNU General Public License for more details.
20
21 You should have received a copy of the GNU General Public License along
22 with this program; if not, write to the Free Software Foundation, Inc.,
23 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24
25 See the full license in the file "LICENSE" in the top level distribution directory
26*************************************************************************************/
27/* END LEGAL */
28
29#include <Grid/GridCore.h>
30
32
33// static data
34
36uint64_t GlobalSharedMemory::MAX_MPI_SHM_BYTES = 1024LL*1024LL*1024LL;
41
42std::vector<void *> GlobalSharedMemory::WorldShmCommBufs;
43#ifndef ACCELERATOR_AWARE_MPI
44void * GlobalSharedMemory::HostCommBuf;
45#endif
46
51
55
58
60{
61 assert(_ShmAlloc);
62 assert(_ShmAllocBytes>0);
63 for(int r=0;r<WorldShmSize;r++){
65 }
66 _ShmAlloc = 0;
68}
69
70// Alloc, free shmem region
72#ifndef ACCELERATOR_AWARE_MPI
73void *SharedMemory::HostBufferMalloc(size_t bytes){
74 void *ptr = (void *)host_heap_top;
75 host_heap_top += bytes;
76 host_heap_bytes+= bytes;
77 if (host_heap_bytes >= host_heap_size) {
78 std::cout<< " HostBufferMalloc exceeded heap size -- try increasing with --shm <MB> flag" <<std::endl;
79 std::cout<< " Parameter specified in units of MB (megabytes) " <<std::endl;
80 std::cout<< " Current alloc is " << (bytes/(1024*1024)) <<"MB"<<std::endl;
81 std::cout<< " Current bytes is " << (host_heap_bytes/(1024*1024)) <<"MB"<<std::endl;
82 std::cout<< " Current heap is " << (host_heap_size/(1024*1024)) <<"MB"<<std::endl;
83 assert(host_heap_bytes<host_heap_size);
84 }
85 return ptr;
86}
87void SharedMemory::HostBufferFreeAll(void) {
88 host_heap_top =(size_t)HostCommBuf;
89 host_heap_bytes=0;
90}
91#endif
93 // bytes = (bytes+sizeof(vRealD))&(~(sizeof(vRealD)-1));// align up bytes
94 void *ptr = (void *)heap_top;
95 heap_top += bytes;
96 heap_bytes+= bytes;
97 if (heap_bytes >= heap_size) {
98 std::cout<< " ShmBufferMalloc exceeded shared heap size -- try increasing with --shm <MB> flag" <<std::endl;
99 std::cout<< " Parameter specified in units of MB (megabytes) " <<std::endl;
100 std::cout<< " Current alloc is " << (bytes/(1024*1024)) <<"MB"<<std::endl;
101 std::cout<< " Current bytes is " << (heap_bytes/(1024*1024)) <<"MB"<<std::endl;
102 std::cout<< " Current heap is " << (heap_size/(1024*1024)) <<"MB"<<std::endl;
103 assert(heap_bytes<heap_size);
104 }
105 //std::cerr << "ShmBufferMalloc "<<std::hex<< ptr<<" - "<<((uint64_t)ptr+bytes)<<std::dec<<std::endl;
106 return ptr;
107}
109 heap_top =(size_t)ShmBufferSelf();
110 heap_bytes=0;
111}
113{
114 //std::cerr << "ShmBufferSelf "<<ShmRank<<" "<<std::hex<< ShmCommBufs[ShmRank] <<std::dec<<std::endl;
115 return ShmCommBufs[ShmRank];
116}
117static inline int divides(int a,int b)
118{
119 return ( b == ( (b/a)*a ) );
120}
122{
124 // Allow user to configure through environment variable
126 char* str = getenv(("GRID_SHM_DIMS_" + std::to_string(ShmDims.size())).c_str());
127 if ( str ) {
128 std::vector<int> IntShmDims;
129 GridCmdOptionIntVector(std::string(str),IntShmDims);
130 assert(IntShmDims.size() == WorldDims.size());
131 long ShmSize = 1;
132 for (int dim=0;dim<WorldDims.size();dim++) {
133 ShmSize *= (ShmDims[dim] = IntShmDims[dim]);
134 assert(divides(ShmDims[dim],WorldDims[dim]));
135 }
136 assert(ShmSize == WorldShmSize);
137 return;
138 }
139
141 // Powers of 2,3,5 only in prime decomposition for now
143 int ndimension = WorldDims.size();
144 ShmDims=Coordinate(ndimension,1);
145
146 std::vector<int> primes({2,3,5});
147
148 int dim = 0;
149 int last_dim = ndimension - 1;
150 int AutoShmSize = 1;
151 while(AutoShmSize != WorldShmSize) {
152 int p;
153 for(p=0;p<primes.size();p++) {
154 int prime=primes[p];
155 if ( divides(prime,WorldDims[dim]/ShmDims[dim])
156 && divides(prime,WorldShmSize/AutoShmSize) ) {
157 AutoShmSize*=prime;
158 ShmDims[dim]*=prime;
159 last_dim = dim;
160 break;
161 }
162 }
163 if (p == primes.size() && last_dim == dim) {
164 std::cerr << "GlobalSharedMemory::GetShmDims failed" << std::endl;
165 exit(EXIT_FAILURE);
166 }
167 dim=(dim+1) %ndimension;
168 }
169}
170
172
AcceleratorVector< int, MaxDims > Coordinate
Definition Coordinate.h:95
void GridCmdOptionIntVector(const std::string &str, VectorInt &vec)
Definition Init.cc:160
#define NAMESPACE_BEGIN(A)
Definition Namespace.h:35
#define NAMESPACE_END(A)
Definition Namespace.h:36
static int divides(int a, int b)
int Grid_MPI_Comm
accelerator_inline size_type size(void) const
Definition Coordinate.h:52
static void SharedMemoryFree(void)
static int WorldShmSize
static Grid_MPI_Comm WorldComm
static uint64_t _ShmAllocBytes
static std::vector< int > WorldShmRanks
static void GetShmDims(const Coordinate &WorldDims, Coordinate &ShmDims)
static int _ShmAlloc
static uint64_t MAX_MPI_SHM_BYTES
static int _ShmSetup
static int HPEhypercube
static int WorldNodes
static int WorldShmRank
static std::vector< void * > WorldShmCommBufs
static Grid_MPI_Comm WorldShmComm
void * ShmBufferSelf(void)
void * ShmBufferMalloc(size_t bytes)
size_t heap_bytes
void ShmBufferFreeAll(void)
std::vector< void * > ShmCommBufs