NAMD
MigrationCUDAKernel.h
Go to the documentation of this file.
1 #ifndef MIGRATIONCUDA_H
2 #define MIGRATIONCUDA_H
3 #include "NamdTypes.h"
4 #include "CudaUtils.h"
5 #include "CudaRecord.h"
6 #include "Lattice.h"
7 #ifdef NODEGROUP_FORCE_REGISTER
8 
9 class MigrationCUDAKernel {
10 private:
11  int* d_patchOffset_temp;
12  int* d_patchOffsetNB_temp;
13  int patchDeviceScan_alloc;
14  char* d_patchDeviceScan_scratch;
15 
16 public:
17  static constexpr int kSortNumThreads = 512;
18  static constexpr int kValuesPerThread = 4;
19  static constexpr int kMaxAtomsPerPatch = kSortNumThreads * kValuesPerThread;
20 
21  static constexpr int kAtomsPerBuffer = 32;
22  static constexpr int kNumSoABuffers = 17;
23 
24  MigrationCUDAKernel();
25  ~MigrationCUDAKernel();
26  void allocateScratch(const int numPatchesHomeAndProxy);
27 
28  void sortAtoms(
29  const int numPatches,
30  const int numAtoms,
31  const CudaLocalRecord* records,
32  const double3* patchMin,
33  const double3* patchMax,
34  const double* pos_x,
35  const double* pos_y,
36  const double* pos_z,
37  int* sortOrder,
38  int* unsortOrder,
39  int* sortIndex,
40  cudaStream_t stream
41  );
42 
43  void copy_AoS_to_SoA(
44  const int numPatches,
45  const bool alchOn,
46  const bool langevinOn,
47  const double dt,
48  const double kbT,
49  const double tempFactor,
50  const CudaLocalRecord* records,
51  const FullAtom* atomdata_AoS,
52  double* recipMass,
53  double* vel_x,
54  double* vel_y,
55  double* vel_z,
56  double* pos_x,
57  double* pos_y,
58  double* pos_z,
59  float* mass,
60  float* charge,
61  int* id,
62  int* vdwType,
63  int* hydrogenGroupSize,
64  int* migrationGroupSize,
65  int* atomFixed,
66  float* rigidBondLength,
67  char3* transform,
68  int* partition,
69  float* langevinParam,
70  float* langScalVelBBK2,
71  float* langScalRandBBK2,
72  cudaStream_t stream
73  );
74 
75  void sortSolventAtoms(
76  const int numPatches,
77  const CudaLocalRecord* records,
78  const FullAtom* atomdata_AoS_in,
79  FullAtom* atomdata_AoS_out,
80  int* sortIndex,
81  cudaStream_t stream
82  );
83 
84  void computeMigrationGroupIndex(
85  const int numPatches,
86  CudaLocalRecord* records,
87  const int* migrationGroupSize,
88  int* migrationGroupIndex,
89  cudaStream_t stream
90  );
91 
92  void transformMigratedPositions(
93  const int numPatches,
94  const CudaLocalRecord* localRecords,
95  const double3* patchCenter,
96  FullAtom* atomdata_AoS,
97  const Lattice lattice,
98  cudaStream_t stream
99  );
100 
101  void transformPositions(
102  const int numPatches,
103  const CudaLocalRecord* localRecords,
104  const double3* patchCenter,
105  FullAtom* atomdata_AoS,
106  const Lattice lattice,
107  const int* hydrogenGroupSize,
108  const int* migrationGroupSize,
109  const int* migrationGroupIndex,
110  double* pos_x,
111  double* pos_y,
112  double* pos_z,
113  cudaStream_t stream
114  );
115 
116  void computeMigrationDestination(
117  const int numPatches,
118  CudaLocalRecord* localRecords,
119  const Lattice lattice,
120  const CudaMInfo* mInfo,
121  const int* patchToDeviceMap,
122  const int* globalToLocalMap,
123  const double3* patchMin,
124  const double3* patchMax,
125  const int* hydrogenGroupSize,
126  const int* migrationGroupSize,
127  const int* migrationGroupIndex,
128  const double* pos_x,
129  const double* pos_y,
130  const double* pos_z,
131  int4* migrationDestination,
132  cudaStream_t stream
133  );
134 
135  void update_AoS(
136  const int numPatches,
137  const CudaLocalRecord* records,
138  FullAtom* atomdata_AoS,
139  const double* vel_x,
140  const double* vel_y,
141  const double* vel_z,
142  const double* pos_x,
143  const double* pos_y,
144  const double* pos_z,
145  cudaStream_t stream
146  );
147 
148  void performLocalMigration(
149  const int numPatches,
150  CudaLocalRecord* records,
151  const FullAtom* atomdata_AoS_in,
152  FullAtom* atomdata_AoS_out,
153  int4* migrationDestination,
154  cudaStream_t stream
155  );
156 
157  void performMigration(
158  const int numPatches,
159  CudaLocalRecord* records,
160  CudaLocalRecord** peer_records,
161  const FullAtom* local_atomdata_AoS,
162  FullAtom** peer_atomdata_AoS,
163  const int* migrationGroupSize,
164  const int* migrationGroupIndex,
165  int4* migrationDestination,
166  cudaStream_t stream
167  );
168 
169  void updateMigrationDestination(
170  const int numAtomsHome,
171  int4* migrationDestination,
172  int** d_peer_sortSoluteIndex,
173  cudaStream_t stream
174  );
175 
176  void copyDataToProxies(
177  const int deviceIndex,
178  const int numPatchesHome,
179  const int numPatchesHomeAndProxy,
180  const CudaLocalRecord* records,
181  int** peer_id,
182  int** peer_vdwType,
183  int** peer_sortOrder,
184  int** peer_unsortOrder,
185  float** peer_charge,
186  int** peer_partition,
187  double3** peer_patchCenter,
188  bool doAlch,
189  cudaStream_t stream
190  );
191 
192  void copyMigrationDestinationToProxies(
193  const int deviceIndex,
194  const int numPatchesHome,
195  const int numPatchesHomeAndProxy,
196  const CudaLocalRecord* records,
197  const CudaPeerRecord* peerRecords,
198  int4** peer_migrationDestination,
199  cudaStream_t stream
200  );
201 
202  void updateLocalRecords(
203  const int numPatchesHome,
204  CudaLocalRecord* records,
205  CudaLocalRecord** peer_records,
206  const CudaPeerRecord* peerRecords,
207  cudaStream_t stream
208  );
209 
210 
211  void updateLocalRecordsOffset(
212  const int numPatchesHomeAndProxy,
213  CudaLocalRecord* records,
214  cudaStream_t stream
215  );
216 
217  void updatePeerRecords(
218  const int numPatchesHomeAndProxy,
219  CudaLocalRecord* records,
220  CudaLocalRecord** peer_records,
221  CudaPeerRecord* peerRecords,
222  cudaStream_t stream
223  );
224 
225 };
226 
227 #endif // NAMD_CUDA
228 #endif // MIGRATIONCUDA_H
static void partition(int *order, const FullAtom *atoms, int begin, int end)
Definition: SortAtoms.C:45