1#define MAX2(x,y) (((x)>(y))?(x):(y))
2#define MAX3(x,y,z) MAX2((x),MAX2((y),(z)))
6__device__
void GetComponentMol(
Atoms* System,
size_t& Mol,
size_t& comp,
size_t Ncomp)
9 for(
size_t ijk = comp; ijk < Ncomp; ijk++)
11 size_t Mol_ijk = System[ijk].size / System[ijk].Molsize;
21__device__
void ScaleCopyPositions(
Atoms* d_a,
Boxsize Box,
size_t comp,
double AScale,
size_t start,
size_t end, double3 COM)
23 for(
size_t atom = 0; atom < d_a[comp].Molsize; atom++)
25 size_t from = start + atom;
27 xyz = d_a[comp].pos[from];
28 double3 posvec = xyz - COM;
29 PBC(posvec, Box.Cell, Box.InverseCell, Box.Cubic);
31 newCOM.x = COM.x * AScale;
32 newCOM.y = COM.y * AScale;
33 newCOM.z = COM.z * AScale;
35 newxyz = newCOM + posvec;
38 size_t to = end + atom;
39 d_a[comp].pos[to] = newxyz;
42__global__
void ScalePositions(
Atoms* d_a,
Boxsize Box,
double AScale,
size_t NComponent,
bool ScaleFramework,
size_t TotalMol,
bool noCharges,
bool* flag)
45 size_t i = blockIdx.x * blockDim.x + threadIdx.x;
48 size_t Mol = i;
size_t comp = 0;
49 if(!ScaleFramework) comp++;
50 size_t start_comp = comp;
53 for(
size_t ijk = start_comp; ijk < NComponent; ijk++)
55 size_t Mol_ijk = d_a[ijk].size / d_a[ijk].Molsize;
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;
67 COM = d_a[comp].pos[start];
68 ScaleCopyPositions(d_a, Box, comp, AScale, start, end, COM);
73 double InverseScale = 1.0 / AScale;
74 for(
size_t j = 0; j < 9; j++)
76 Box.Cell[j] *= AScale; Box.InverseCell[j] *= InverseScale;
78 Box.Volume *= pow(AScale, 3);
84 double Alpha = Box.Alpha;
85 double tol1 = Box.tol1;
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);
97__global__
void CopyScaledPositions(
Atoms* d_a,
size_t NComponent,
bool ScaleFramework,
size_t TotalMol)
99 size_t i = blockIdx.x * blockDim.x + threadIdx.x;
102 size_t Mol = i;
size_t comp = 0;
103 if(!ScaleFramework) comp++;
104 size_t start_comp = comp;
107 for(
size_t ijk = start_comp; ijk < NComponent; ijk++)
109 size_t Mol_ijk = d_a[ijk].size / d_a[ijk].Molsize;
117 size_t HalfAllocateSize = d_a[comp].Allocate_size / 2;
118 size_t start = Mol * d_a[comp].Molsize;
119 size_t end = start + HalfAllocateSize;
121 for(
size_t atom = 0; atom < d_a[comp].Molsize; atom++)
123 size_t from = end + atom;
124 size_t to = start + atom;
125 d_a[comp].pos[to] = d_a[comp].pos[from];
130__global__
void Revert_Boxsize(
Boxsize Box,
double AScale,
bool noCharges)
132 size_t i = blockIdx.x * blockDim.x + threadIdx.x;
135 double InverseScale = 1.0 / AScale;
136 for(
size_t j = 0; j < 9; j++)
138 Box.Cell[j] *= InverseScale; Box.InverseCell[j] *= AScale;
140 Box.Volume *= pow(InverseScale, 3);
143 double Alpha = Box.Alpha;
144 double tol1 = Box.tol1;
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);
156static inline double Get_TotalNumberOfMolecule_In_Box(
Components& SystemComponents)
158 double NMol =
static_cast<double>(SystemComponents.TotalNumberOfMolecules - SystemComponents.NumberOfFrameworks);
160 for(
size_t comp = 0; comp < SystemComponents.Total_Components; comp++)
162 if(SystemComponents.hasfractionalMolecule[comp])
170void NVTGibbsMove(std::vector<Components>& SystemComponents,
Simulations*& Sims,
ForceField FF, std::vector<SystemEnergies>& Energy,
Gibbs& GibbsStatistics)
172 size_t NBox = SystemComponents.size();
173 size_t SelectedBox = 0;
175 if(Get_Uniform_Random() > 0.5)
181 GibbsStatistics.GibbsBoxStats.x += 1.0;
182 double TotalV = Sims[SelectedBox].Box.Volume + Sims[OtherBox].Box.Volume;
184 double OldVA = Sims[SelectedBox].Box.Volume;
double OldVB = Sims[OtherBox].Box.Volume;
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;
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;
203 bool Overlap =
false;
207 for(
size_t sim = 0; sim < NBox; sim++)
209 if(Overlap)
continue;
210 totMol = SystemComponents[sim].TotalNumberOfMolecules - SystemComponents[sim].NumberOfFrameworks;
211 for(
size_t comp = 1; comp < SystemComponents[sim].Total_Components; comp++)
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!!!");
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);
223 size_t Host_threads = 0;
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);
230 cudaMemcpy(SystemComponents[sim].flag, Sims[sim].device_flag,
sizeof(
bool), cudaMemcpyDeviceToHost);
231 if(SystemComponents[sim].flag[0])
234 printf(
"There is OVERLAP IN GIBBS VOLUME MOVE in Box[%zu]!\n", sim);
242 bool UseOffSet =
true;
243 NewE[sim] += Ewald_TotalEnergy(Sims[sim], SystemComponents[sim], UseOffSet);
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;
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;
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))));
264 if(Get_Uniform_Random() < Pacc) Accept =
true;
269 GibbsStatistics.GibbsBoxStats.y += 1.0;
271 for(
size_t sim = 0; sim < NBox; sim++)
273 totMol = SystemComponents[sim].TotalNumberOfMolecules - SystemComponents[sim].NumberOfFrameworks;
274 for(
size_t comp = 1; comp < SystemComponents[sim].Total_Components; comp++)
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!!!");
279 double CurrentEnergy = Energy[sim].running_energy + Energy[sim].InitialEnergy;
281 SystemComponents[sim].deltaE += DeltaE[sim];
283 Setup_threadblock(totMol, &Nblock, &Nthread);
285 CopyScaledPositions<<<Nblock, Nthread>>>(Sims[sim].d_a, SystemComponents[sim].Total_Components, ScaleFramework, totMol);
288 Update_Ewald_Vector(Sims[sim].Box,
false, SystemComponents[sim], 0);
293 for(
size_t sim = 0; sim < NBox; sim++)
296 Revert_Boxsize<<<1,1>>>(Sims[sim].Box, ScaleAB[sim], FF.noCharges);
301static inline void Update_Max_GibbsVolume(
Gibbs& GibbsStatistics)
303 if(GibbsStatistics.GibbsBoxStats.x > 0)
305 double ratio =
static_cast<double>(GibbsStatistics.GibbsBoxStats.x) /
static_cast<double>(GibbsStatistics.GibbsBoxStats.y);
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;
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