My Project
Loading...
Searching...
No Matches
mc_box.h
1#define MAX2(x,y) (((x)>(y))?(x):(y)) // the maximum of two numbers
2#define MAX3(x,y,z) MAX2((x),MAX2((y),(z)))
3
4void NVTGibbsMove(std::vector<Components>& SystemComponents, Simulations*& Sims, ForceField FF);
5
6__device__ void GetComponentMol(Atoms* System, size_t& Mol, size_t& comp, size_t Ncomp)
7{
8 size_t total = 0;
9 for(size_t ijk = comp; ijk < Ncomp; ijk++)
10 {
11 size_t Mol_ijk = System[ijk].size / System[ijk].Molsize;
12 total += Mol_ijk;
13 if(Mol >= total)
14 {
15 comp++;
16 Mol -= Mol_ijk;
17 }
18 }
19}
20
21__device__ void ScaleCopyPositions(Atoms* d_a, Boxsize Box, size_t comp, double AScale, size_t start, size_t end, double3 COM)
22{
23 for(size_t atom = 0; atom < d_a[comp].Molsize; atom++)
24 {
25 size_t from = start + atom;
26 double3 xyz;
27 xyz = d_a[comp].pos[from];
28 double3 posvec = xyz - COM;
29 PBC(posvec, Box.Cell, Box.InverseCell, Box.Cubic);
30 double3 newCOM;
31 newCOM.x = COM.x * AScale;
32 newCOM.y = COM.y * AScale;
33 newCOM.z = COM.z * AScale;
34 double3 newxyz;
35 newxyz = newCOM + posvec;
36 //Copy data to new location//
37 //Copy to the later half of the xyz arrays in d_a//
38 size_t to = end + atom;
39 d_a[comp].pos[to] = newxyz;
40 }
41}
42__global__ void ScalePositions(Atoms* d_a, Boxsize Box, double AScale, size_t NComponent, bool ScaleFramework, size_t TotalMol, bool noCharges, bool* flag)
43{
44 //Each thread processes an adsorbate molecule//
45 size_t i = blockIdx.x * blockDim.x + threadIdx.x;
46 if(i < TotalMol)
47 {
48 size_t Mol = i; size_t comp = 0;
49 if(!ScaleFramework) comp++;
50 size_t start_comp = comp;
51 //GetComponentMol(d_a, Mol, comp, NComponent);
52 size_t total = 0;
53 for(size_t ijk = start_comp; ijk < NComponent; ijk++)
54 {
55 size_t Mol_ijk = d_a[ijk].size / d_a[ijk].Molsize;
56 total += Mol_ijk;
57 if(Mol >= total)
58 {
59 comp++;
60 Mol -= Mol_ijk;
61 }
62 }
63 size_t HalfAllocateSize = d_a[comp].Allocate_size / 2;
64 size_t start = Mol * d_a[comp].Molsize;
65 size_t end = start + HalfAllocateSize;
66 double3 COM; //Zhao's note: one can decide whether to use the first bead position or the center of mass//
67 COM = d_a[comp].pos[start];
68 ScaleCopyPositions(d_a, Box, comp, AScale, start, end, COM);
69 }
70 //Scale Boxsize//
71 if(i == 0)
72 {
73 double InverseScale = 1.0 / AScale;
74 for(size_t j = 0; j < 9; j++)
75 {
76 Box.Cell[j] *= AScale; Box.InverseCell[j] *= InverseScale;
77 }
78 Box.Volume *= pow(AScale, 3);
79 //Update kmax and reciprocal cutoff//
80 //Only Alpha and tol1 are needed to pass here//
81 //see "read_Ewald_Parameters_from_input" function in read_data.cpp//
82 if(!noCharges)
83 {
84 double Alpha = Box.Alpha;
85 double tol1 = Box.tol1;
86 //Zhao's note: See InitializeEwald function in RASPA-2.0 //
87 // 1.0 / Pi = 0.31830988618 //
88 Box.kmax.x = std::round(0.25 + Box.Cell[0] * Alpha * tol1 * 0.31830988618);
89 Box.kmax.y = std::round(0.25 + Box.Cell[4] * Alpha * tol1 * 0.31830988618);
90 Box.kmax.z = std::round(0.25 + Box.Cell[8] * Alpha * tol1 * 0.31830988618);
91 Box.ReciprocalCutOff = pow(1.05*static_cast<double>(MAX3(Box.kmax.x, Box.kmax.y, Box.kmax.z)), 2);
92 }
93 flag[0] = false;
94 }
95}
96
97__global__ void CopyScaledPositions(Atoms* d_a, size_t NComponent, bool ScaleFramework, size_t TotalMol)
98{
99 size_t i = blockIdx.x * blockDim.x + threadIdx.x;
100 if(i < TotalMol)
101 {
102 size_t Mol = i; size_t comp = 0;
103 if(!ScaleFramework) comp++;
104 size_t start_comp = comp;
105 //GetComponentMol(d_a, Mol, comp, NComponent);
106 size_t total = 0;
107 for(size_t ijk = start_comp; ijk < NComponent; ijk++)
108 {
109 size_t Mol_ijk = d_a[ijk].size / d_a[ijk].Molsize;
110 total += Mol_ijk;
111 if(Mol >= total)
112 {
113 comp++;
114 Mol -= Mol_ijk;
115 }
116 }
117 size_t HalfAllocateSize = d_a[comp].Allocate_size / 2;
118 size_t start = Mol * d_a[comp].Molsize; //position of the original xyz data before the Gibbs move//
119 size_t end = start + HalfAllocateSize; //position of the scaled xyz data after the Gibbs move//
120 //Copy the xyz data from the scaled region (later half) to the original region (first half)//
121 for(size_t atom = 0; atom < d_a[comp].Molsize; atom++)
122 {
123 size_t from = end + atom;
124 size_t to = start + atom;
125 d_a[comp].pos[to] = d_a[comp].pos[from];
126 }
127 }
128}
129
130__global__ void Revert_Boxsize(Boxsize Box, double AScale, bool noCharges)
131{
132 size_t i = blockIdx.x * blockDim.x + threadIdx.x;
133 if(i == 0)
134 {
135 double InverseScale = 1.0 / AScale;
136 for(size_t j = 0; j < 9; j++)
137 {
138 Box.Cell[j] *= InverseScale; Box.InverseCell[j] *= AScale;
139 }
140 Box.Volume *= pow(InverseScale, 3);
141 if(!noCharges)
142 {
143 double Alpha = Box.Alpha;
144 double tol1 = Box.tol1;
145 //Zhao's note: See InitializeEwald function in RASPA-2.0 //
146 // 1.0 / Pi = 0.31830988618 //
147 Box.kmax.x = std::round(0.25 + Box.Cell[0] * Alpha * tol1 * 0.31830988618);
148 Box.kmax.y = std::round(0.25 + Box.Cell[4] * Alpha * tol1 * 0.31830988618);
149 Box.kmax.z = std::round(0.25 + Box.Cell[8] * Alpha * tol1 * 0.31830988618);
150 Box.ReciprocalCutOff = pow(1.05*static_cast<double>(MAX3(Box.kmax.x, Box.kmax.y, Box.kmax.z)), 2);
151 }
152 }
153}
154
155//Get the total number of molecules (excluding framework and fractional molecule) in a box//
156static inline double Get_TotalNumberOfMolecule_In_Box(Components& SystemComponents)
157{
158 double NMol = static_cast<double>(SystemComponents.TotalNumberOfMolecules - SystemComponents.NumberOfFrameworks);
159 //Minus the fractional molecule//
160 for(size_t comp = 0; comp < SystemComponents.Total_Components; comp++)
161 {
162 if(SystemComponents.hasfractionalMolecule[comp])
163 {
164 NMol-=1.0;
165 }
166 }
167 return NMol;
168}
169
170void NVTGibbsMove(std::vector<Components>& SystemComponents, Simulations*& Sims, ForceField FF, std::vector<SystemEnergies>& Energy, Gibbs& GibbsStatistics)
171{
172 size_t NBox = SystemComponents.size();
173 size_t SelectedBox = 0;
174 size_t OtherBox = 1;
175 if(Get_Uniform_Random() > 0.5)
176 {
177 SelectedBox = 1;
178 OtherBox = 0;
179 }
180 //printf("Performing NVT Gibbs Move! on Box[%zu]\n", SelectedBox);
181 GibbsStatistics.GibbsBoxStats.x += 1.0;
182 double TotalV = Sims[SelectedBox].Box.Volume + Sims[OtherBox].Box.Volume;
183
184 double OldVA = Sims[SelectedBox].Box.Volume; double OldVB = Sims[OtherBox].Box.Volume;
185
186 double MaxGibbsVolumeChange = GibbsStatistics.MaxGibbsBoxChange;
187 double expdV = std::exp(std::log(Sims[SelectedBox].Box.Volume / Sims[OtherBox].Box.Volume) + MaxGibbsVolumeChange * 2.0 * (Get_Uniform_Random() - 0.5));
188 double newVA = expdV * TotalV / (1.0 + expdV);
189 double newVB = TotalV - newVA;
190 //printf("TotalV: %.5f, VA/VB: %.5f/%.5f, expdV: %.5f, newVA/newVB: %.5f/%.5f\n", TotalV, Sims[SelectedBox].Box.Volume, Sims[OtherBox].Box.Volume, expdV, newVA, newVB);
191 double ScaleAB[2] = {0.0, 0.0};
192 ScaleAB[SelectedBox] = std::cbrt(newVA / Sims[SelectedBox].Box.Volume);
193 ScaleAB[OtherBox] = std::cbrt(newVB / Sims[OtherBox].Box.Volume);
194 bool ScaleFramework = false;
195 size_t Nblock = 0; size_t Nthread = 0; size_t totMol = 0;
196 //Check if Allocate_size is greater than or equal to twice of the current size//
197 //Then scale the boxes, calculate the new energy//
198
199 MoveEnergy CurrentE[2];
200 MoveEnergy NewE[2];
201 MoveEnergy DeltaE[2];
202
203 bool Overlap = false;
204
205 if(!ScaleFramework)
206 {
207 for(size_t sim = 0; sim < NBox; sim++)
208 {
209 if(Overlap) continue;
210 totMol = SystemComponents[sim].TotalNumberOfMolecules - SystemComponents[sim].NumberOfFrameworks;
211 for(size_t comp = 1; comp < SystemComponents[sim].Total_Components; comp++)
212 {
213 size_t TotSize = SystemComponents[sim].Moleculesize[comp] * SystemComponents[sim].NumberOfMolecule_for_Component[comp];
214 if(TotSize * 2 > SystemComponents[sim].Allocate_size[comp]) throw std::runtime_error("Allocate More space for adsorbates on the GPU!!!");
215 //printf("Box[%zu], TotalMolecule for component [%zu] is %zu\n", sim, comp, totMol);
216 }
217 Setup_threadblock(totMol, &Nblock, &Nthread);
218 ScalePositions<<<Nblock, Nthread>>>(Sims[sim].d_a, Sims[sim].Box, ScaleAB[sim], SystemComponents[sim].Total_Components, ScaleFramework, totMol, FF.noCharges, Sims[sim].device_flag);
219
221 // TOTAL VDW + REAL //
223 size_t Host_threads = 0; //In this move, there should be NO Framework Atom (Framework (component 0) should be an empty box)//
224 size_t Guest_threads = totMol * (totMol - 1) / 2;
225 bool ConsiderHostHost = false;
226 bool UseOffset = true;
227 size_t NFrameworkAtomsPerThread = 1;
228 NewE[sim] = Total_VDW_Coulomb_Energy(Sims[sim], FF, totMol, Host_threads, Guest_threads, NFrameworkAtomsPerThread, ConsiderHostHost, UseOffset);
229 //Check for Overlaps//
230 cudaMemcpy(SystemComponents[sim].flag, Sims[sim].device_flag, sizeof(bool), cudaMemcpyDeviceToHost);
231 if(SystemComponents[sim].flag[0])
232 {
233 Overlap = true;
234 printf("There is OVERLAP IN GIBBS VOLUME MOVE in Box[%zu]!\n", sim);
235 }
237 // TOTAL EWALD CALCULATION //
239 //Calculate new alpha and kmax//
240 if(!FF.noCharges)
241 {
242 bool UseOffSet = true; //Use the second half of the allocated space for xyz
243 NewE[sim] += Ewald_TotalEnergy(Sims[sim], SystemComponents[sim], UseOffSet);
244 }
245 }
246 }
247 //If the Gibbs Volume Move is Accepted, for the test, assume it is always accepted//
248 bool Accept = false;
249 if(!Overlap)
250 {
251 double NMolA= Get_TotalNumberOfMolecule_In_Box(SystemComponents[SelectedBox]);
252 double NMolB= Get_TotalNumberOfMolecule_In_Box(SystemComponents[OtherBox]);
253 CurrentE[SelectedBox] = SystemComponents[SelectedBox].CreateMol_Energy + SystemComponents[SelectedBox].deltaE;
254 CurrentE[OtherBox] = SystemComponents[OtherBox].CreateMol_Energy + SystemComponents[OtherBox].deltaE;
255
256 DeltaE[SelectedBox] = NewE[SelectedBox] - CurrentE[SelectedBox];
257 DeltaE[OtherBox] = NewE[OtherBox] - CurrentE[OtherBox];
258 double VolumeRatioA= Sims[SelectedBox].Box.Volume / OldVA;
259 double VolumeRatioB= Sims[OtherBox].Box.Volume / OldVB;
260
261 //This assumes that the two boxes share the same temperature, it might not be true//
262 double Pacc = std::exp(-SystemComponents[SelectedBox].Beta*((DeltaE[0].total()+DeltaE[1].total())+((NMolA+1.0)*std::log(VolumeRatioA))+((NMolB+1.0)*std::log(VolumeRatioB))));
263 //printf("Gibbs Box Move, Pacc: %.5f\n", Pacc);
264 if(Get_Uniform_Random() < Pacc) Accept = true;
265 }
266
267 if(Accept)
268 {
269 GibbsStatistics.GibbsBoxStats.y += 1.0;
270 //Update Energy and positions//
271 for(size_t sim = 0; sim < NBox; sim++)
272 {
273 totMol = SystemComponents[sim].TotalNumberOfMolecules - SystemComponents[sim].NumberOfFrameworks;
274 for(size_t comp = 1; comp < SystemComponents[sim].Total_Components; comp++)
275 {
276 size_t TotSize = SystemComponents[sim].Moleculesize[comp] * SystemComponents[sim].NumberOfMolecule_for_Component[comp];
277 if(TotSize * 2 > SystemComponents[sim].Allocate_size[comp]) throw std::runtime_error("Allocate More space for adsorbates on the GPU!!!");
278 }
279 double CurrentEnergy = Energy[sim].running_energy + Energy[sim].InitialEnergy;
280 //Zhao's note: here the energies are only considering VDW + Real//
281 SystemComponents[sim].deltaE += DeltaE[sim];
282 //printf("OldEnergy: %.5f, NewEnergy: %.5f, DeltaE: %.5f\n", CurrentEnergy, NewEnergy, DeltaE);
283 Setup_threadblock(totMol, &Nblock, &Nthread);
284 //Copy xyz data from new to old, also update box lengths//
285 CopyScaledPositions<<<Nblock, Nthread>>>(Sims[sim].d_a, SystemComponents[sim].Total_Components, ScaleFramework, totMol);
286 //Update Eik if accepted from tempEik to StoredEik, BUG (adsorbate/framework species all needs to be updated)!!!//
287 if(!FF.noCharges)
288 Update_Ewald_Vector(Sims[sim].Box, false, SystemComponents[sim], 0);
289 }
290 }
291 else
292 {
293 for(size_t sim = 0; sim < NBox; sim++)
294 {
295 //Revert the Boxsize, if Charges, update kmax and Reciprocal Cutoff//
296 Revert_Boxsize<<<1,1>>>(Sims[sim].Box, ScaleAB[sim], FF.noCharges);
297 }
298 }
299}
300
301static inline void Update_Max_GibbsVolume(Gibbs& GibbsStatistics)
302{
303 if(GibbsStatistics.GibbsBoxStats.x > 0)
304 {
305 double ratio = static_cast<double>(GibbsStatistics.GibbsBoxStats.x) / static_cast<double>(GibbsStatistics.GibbsBoxStats.y);
306 double vandr = 1.05;
307 if(ratio < 0.5) vandr = 0.95;
308 GibbsStatistics.MaxGibbsBoxChange*=vandr;
309 if(GibbsStatistics.MaxGibbsBoxChange<0.0005)
310 GibbsStatistics.MaxGibbsBoxChange=0.0005;
311 if(GibbsStatistics.MaxGibbsBoxChange>0.5)
312 GibbsStatistics.MaxGibbsBoxChange=0.5;
313 }
314}
Definition data_struct.h:746
Definition data_struct.h:819
Definition data_struct.h:843
Definition data_struct.h:794
Definition data_struct.h:61
Definition data_struct.h:615
Definition data_struct.h:1027