Swarm-NG  1.1
bppt.hpp
Go to the documentation of this file.
1 /*************************************************************************
2  * Copyright (C) 2011 by Saleh Dindar and the Swarm-NG Development Team *
3  * *
4  * This program is free software; you can redistribute it and/or modify *
5  * it under the terms of the GNU General Public License as published by *
6  * the Free Software Foundation; either version 3 of the License. *
7  * *
8  * This program is distributed in the hope that it will be useful, *
9  * but WITHOUT ANY WARRANTY; without even the implied warranty of *
10  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
11  * GNU General Public License for more details. *
12  * *
13  * You should have received a copy of the GNU General Public License *
14  * along with this program; if not, write to the *
15  * Free Software Foundation, Inc., *
16  * 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA. *
17  ************************************************************************/
18 
27 #pragma once
28 
29 #include "../integrator.hpp"
30 #include "../log/log.hpp"
31 #include "../plugin.hpp"
32 
33 #include "helpers.hpp"
34 #include "utilities.hpp"
35 #include "device_settings.hpp"
36 
37 
38 namespace swarm {
39 namespace gpu {
40 
68 namespace bppt {
69 
70 
74 GPUAPI int sysid(){
75  return ((blockIdx.z * gridDim.y + blockIdx.y) * gridDim.x + blockIdx.x) * blockDim.x + threadIdx.x;
76 }
77 
81 GPUAPI int sysid_in_block(){
82  return threadIdx.x;
83 }
84 
88 GPUAPI int thread_in_system() {
89  return threadIdx.y;
90 }
91 
92 
96 GPUAPI int system_per_block_gpu() {
97  return blockDim.x;
98  };
99 
103 GPUAPI int thread_component_idx(int nbod) {
104  return thread_in_system() / nbod;
105 }
106 
110 GPUAPI int thread_body_idx(int nbod) {
111  return thread_in_system() % nbod;
112 }
113 
114 
121 template< class Impl, class T>
122 GPUAPI void * system_shared_data_pointer(Impl* integ, T compile_time_param) {
123  extern __shared__ char shared_mem[];
124  int b = sysid_in_block() / SHMEM_CHUNK_SIZE ;
125  int i = sysid_in_block() % SHMEM_CHUNK_SIZE ;
126  int idx = i * sizeof(double)
127  + b * SHMEM_CHUNK_SIZE
128  * Impl::shmem_per_system(compile_time_param);
129  return &shared_mem[idx];
130 }
131 
132 
133 
140 class integrator : public gpu::integrator {
141  typedef gpu::integrator Base;
142  protected:
146 
148  public:
149  integrator(const config &cfg) : Base(cfg) {
150  int spb = cfg.optional("system_per_block",0);
151  _override_system_per_block =(spb % SHMEM_CHUNK_SIZE == 0) ? spb : 0;
152  }
153 
155 
162  const int& override_system_per_block()const{
164  }
165 
177  template<class T>
178  static GENERIC int thread_per_system(T compile_time_param){
179  const int nbod = T::n;
180  const int body_comp = nbod * 3;
181  const int pair_count = nbod * (nbod - 1) / 2;
182  return (body_comp>pair_count) ? body_comp : pair_count;
183  }
184 
185 
194  template<class T>
195  static GENERIC int shmem_per_system(T compile_time_param){
196  const int nbod = T::n;
197  const int pair_count = nbod * (nbod - 1) / 2;
198  return pair_count * 3 * 2 * sizeof(double);
199  }
200 
201 
202 };
203 
204 
205 }
206 }
207 }
208