My Project
Loading...
Searching...
No Matches
mc_swap_moves.h
1#include "mc_widom.h"
2#include "mc_swap_utilities.h"
3#include "lambda.h"
4#include "mc_cbcfc.h"
5
6__global__ void StoreNewLocation_Reinsertion(Atoms Mol, Atoms NewMol, double3* temp, size_t SelectedTrial, size_t Moleculesize)
7{
8 if(Moleculesize == 1) //Only first bead is inserted, first bead data is stored in NewMol
9 {
10 temp[0] = NewMol.pos[SelectedTrial];
11 }
12 else //Multiple beads: first bead + trial orientations
13 {
14 //Update the first bead, first bead data stored in position 0 of Mol //
15 temp[0] = Mol.pos[0];
16
17 size_t chainsize = Moleculesize - 1; // FOr trial orientations //
18 for(size_t i = 0; i < chainsize; i++) //Update the selected orientations//
19 {
20 size_t selectsize = SelectedTrial*chainsize+i;
21 temp[i+1] = NewMol.pos[selectsize];
22 }
23 }
24 /*
25 for(size_t i = 0; i < Moleculesize; i++)
26 printf("i: %lu, xyz: %.5f %.5f %.5f\n", i, temp[i].x, temp[i].y, temp[i].z);
27 */
28}
29
30__global__ void Update_Reinsertion_data(Atoms* d_a, double3* temp, size_t SelectedComponent, size_t UpdateLocation)
31{
32 size_t i = blockIdx.x * blockDim.x + threadIdx.x;
33 size_t realLocation = UpdateLocation + i;
34 d_a[SelectedComponent].pos[realLocation] = temp[i];
35}
36
37static inline MoveEnergy Reinsertion(Components& SystemComponents, Simulations& Sims, ForceField& FF, RandomNumber& Random, WidomStruct& Widom, size_t SelectedMolInComponent, size_t SelectedComponent)
38{
39 //Get Number of Molecules for this component (For updating TMMC)//
40 double NMol = SystemComponents.NumberOfMolecule_for_Component[SelectedComponent];
41 if(SystemComponents.hasfractionalMolecule[SelectedComponent]) NMol--;
42
43 SystemComponents.Moves[SelectedComponent].ReinsertionTotal ++;
44 bool SuccessConstruction = false;
45 size_t SelectedTrial = 0;
46 MoveEnergy energy; MoveEnergy old_energy; double StoredR = 0.0;
47
49 // INSERTION //
51 int CBMCType = REINSERTION_INSERTION; //Reinsertion-Insertion//
52 double2 newScale = SystemComponents.Lambda[SelectedComponent].SET_SCALE(1.0); //Zhao's note: not used in reinsertion, just set to 1.0//
53 double Rosenbluth=Widom_Move_FirstBead_PARTIAL(SystemComponents, Sims, FF, Random, Widom, SelectedMolInComponent, SelectedComponent, CBMCType, StoredR, &SelectedTrial, &SuccessConstruction, &energy, newScale); //Not reinsertion, not Retrace//
54
55 if(Rosenbluth <= 1e-150) SuccessConstruction = false; //Zhao's note: added this protection bc of weird error when testing GibbsParticleXfer
56
57 if(!SuccessConstruction)
58 {
59 SystemComponents.Tmmc[SelectedComponent].Update(1.0, NMol, REINSERTION);
60 energy.zero();
61 return energy;
62 }
63
64 //DEBUG//
65 /*
66 if(SystemComponents.CURRENTCYCLE == 28)
67 {printf("REINSERTION MOVE (comp: %zu, Mol: %zu) FIRST BEAD INSERTION ENERGY: ", SelectedComponent, SelectedMolInComponent); energy.print();
68 printf("Rosen: %.5f\n", Rosenbluth);
69 }
70 */
71
72 if(SystemComponents.Moleculesize[SelectedComponent] > 1)
73 {
74 size_t SelectedFirstBeadTrial = SelectedTrial;
75 MoveEnergy temp_energy = energy;
76 Rosenbluth*=Widom_Move_Chain_PARTIAL(SystemComponents, Sims, FF, Random, Widom, SelectedMolInComponent, SelectedComponent, CBMCType, &SelectedTrial, &SuccessConstruction, &energy, SelectedFirstBeadTrial, newScale); //True for doing insertion for reinsertion, different in MoleculeID//
77 if(Rosenbluth <= 1e-150) SuccessConstruction = false; //Zhao's note: added this protection bc of weird error when testing GibbsParticleXfer
78 if(!SuccessConstruction)
79 {
80 SystemComponents.Tmmc[SelectedComponent].Update(1.0, NMol, REINSERTION);
81 energy.zero();
82 return energy;
83 }
84 energy += temp_energy;
85 }
86
87 /*
88 if(SystemComponents.CURRENTCYCLE == 28)
89 {printf("REINSERTION MOVE, INSERTION ENERGY: "); energy.print();
90 printf("Rosen: %.5f\n", Rosenbluth);
91 }
92 */
93
94 //Store The New Locations//
95 double3* temp;
96 cudaMalloc(&temp, sizeof(double3) * SystemComponents.Moleculesize[SelectedComponent]);
97 StoreNewLocation_Reinsertion<<<1,1>>>(Sims.Old, Sims.New, temp, SelectedTrial, SystemComponents.Moleculesize[SelectedComponent]);
99 // RETRACE //
101 CBMCType = REINSERTION_RETRACE; //Reinsertion-Retrace//
102 double Old_Rosen=Widom_Move_FirstBead_PARTIAL(SystemComponents, Sims, FF, Random, Widom, SelectedMolInComponent, SelectedComponent, CBMCType, StoredR, &SelectedTrial, &SuccessConstruction, &old_energy, newScale);
103
104 /*
105 if(SystemComponents.CURRENTCYCLE == 28)
106 {printf("REINSERTION MOVE (comp: %zu, Mol: %zu) FIRST BEAD DELETION ENERGY: ", SelectedComponent, SelectedMolInComponent); old_energy.print();
107 printf("Rosen: %.5f\n", Old_Rosen);
108 }
109 */
110
111
112 if(SystemComponents.Moleculesize[SelectedComponent] > 1)
113 {
114 size_t SelectedFirstBeadTrial = SelectedTrial;
115 MoveEnergy temp_energy = old_energy;
116 Old_Rosen*=Widom_Move_Chain_PARTIAL(SystemComponents, Sims, FF, Random, Widom, SelectedMolInComponent, SelectedComponent, CBMCType, &SelectedTrial, &SuccessConstruction, &old_energy, SelectedFirstBeadTrial, newScale);
117 old_energy += temp_energy;
118 }
119
120 /*
121 if(SystemComponents.CURRENTCYCLE == 28)
122 {printf("REINSERTION MOVE, DELETION ENERGY: "); old_energy.print();
123 printf("Rosen: %.5f\n", Old_Rosen);
124 }
125 */
126
127 energy -= old_energy;
128
129 //Calculate Ewald//
130 size_t UpdateLocation = SystemComponents.Moleculesize[SelectedComponent] * SelectedMolInComponent;
131
132 bool EwaldPerformed = false;
133 if(!FF.noCharges && SystemComponents.hasPartialCharge[SelectedComponent])
134 {
135 double2 EwaldE = GPU_EwaldDifference_Reinsertion(Sims.Box, Sims.d_a, Sims.Old, temp, FF, Sims.Blocksum, SystemComponents, SelectedComponent, UpdateLocation);
136
137 energy.GGEwaldE = EwaldE.x;
138 energy.HGEwaldE = EwaldE.y;
139 Rosenbluth *= std::exp(-SystemComponents.Beta * (EwaldE.x + EwaldE.y));
140 EwaldPerformed = true;
141 }
142 //Calculate DNN//
143 //Put it after Ewald summation, the required positions are already in place (done by the preparation parts of Ewald Summation)//
144 if(SystemComponents.UseDNNforHostGuest)
145 {
146 if(!EwaldPerformed) Prepare_DNN_InitialPositions_Reinsertion(Sims.d_a, Sims.Old, temp, SystemComponents, SelectedComponent, UpdateLocation);
147 energy.DNN_E = DNN_Prediction_Reinsertion(SystemComponents, Sims, SelectedComponent, temp);
148 //printf("DNN Delta Reinsertion: %.5f\n", energy.DNN_E);
149 //Correction of DNN - HostGuest energy to the Rosenbluth weight//
150 double correction = energy.DNN_Correction();
151 if(fabs(correction) > SystemComponents.DNNDrift) //If there is a huge drift in the energy correction between DNN and Classical HostGuest//
152 {
153 //printf("REINSERTION: Bad Prediction, reject the move!!!\n");
154 SystemComponents.ReinsertionDNNReject++;
155 //WriteOutliers(SystemComponents, Sims, REINSERTION_NEW, energy, correction);
156 //WriteOutliers(SystemComponents, Sims, REINSERTION_OLD, energy, correction);
157 energy.zero();
158 return energy;
159 }
160 SystemComponents.ReinsertionDNNDrift += fabs(correction);
161 Rosenbluth *= std::exp(-SystemComponents.Beta * correction);
162 }
163
164 //Determine whether to accept or reject the insertion
165 double RANDOM = Get_Uniform_Random();
166 //printf("RANDOM: %.5f, Rosenbluth / Old_Rosen: %.5f\n", RANDOM, Rosenbluth / Old_Rosen);
167 if(RANDOM < Rosenbluth / Old_Rosen)
168 { // accept the move
169 SystemComponents.Moves[SelectedComponent].ReinsertionAccepted ++;
170 //size_t UpdateLocation = SystemComponents.Moleculesize[SelectedComponent] * SelectedMolInComponent;
171 Update_Reinsertion_data<<<1,SystemComponents.Moleculesize[SelectedComponent]>>>(Sims.d_a, temp, SelectedComponent, UpdateLocation); checkCUDAError("error Updating Reinsertion data");
172 cudaFree(temp);
173 if(!FF.noCharges && SystemComponents.hasPartialCharge[SelectedComponent])
174 Update_Ewald_Vector(Sims.Box, false, SystemComponents, SelectedComponent);
175 SystemComponents.Tmmc[SelectedComponent].Update(1.0, NMol, REINSERTION); //Update for TMMC, since Macrostate not changed, just add 1.//
176 //energy.print();
177 return energy;
178 }
179 else
180 {
181 cudaFree(temp);
182 SystemComponents.Tmmc[SelectedComponent].Update(1.0, NMol, REINSERTION); //Update for TMMC, since Macrostate not changed, just add 1.//
183 energy.zero();
184 return energy;
185 }
186}
187
188//Zhao's note: added feature for creating fractional molecules//
189static inline MoveEnergy CreateMolecule(Components& SystemComponents, Simulations& Sims, ForceField& FF, RandomNumber& Random, WidomStruct& Widom, size_t SelectedMolInComponent, size_t SelectedComponent, double2 newScale)
190{
191 bool SuccessConstruction = false;
192 double Rosenbluth = 0.0;
193 size_t SelectedTrial = 0;
194 double preFactor = 0.0;
195
196 //Zhao's note: For creating the fractional molecule, there is no previous step, so set the flag to false//
197 MoveEnergy energy = Insertion_Body(SystemComponents, Sims, FF, Random, Widom, SelectedMolInComponent, SelectedComponent, Rosenbluth, SuccessConstruction, SelectedTrial, preFactor, false, newScale);
198 if(!SuccessConstruction)
199 {
200 energy.zero();
201 return energy;
202 }
203 double IdealRosen = SystemComponents.IdealRosenbluthWeight[SelectedComponent];
204 double RANDOM = 1e-100;
205 if(RANDOM < preFactor * Rosenbluth / IdealRosen)
206 { // accept the move
207 size_t UpdateLocation = SystemComponents.Moleculesize[SelectedComponent] * SystemComponents.NumberOfMolecule_for_Component[SelectedComponent];
208 //Zhao's note: here needs more consideration: need to update after implementing polyatomic molecule
209 Update_insertion_data<<<1,1>>>(Sims.d_a, Sims.Old, Sims.New, SelectedTrial, SelectedComponent, UpdateLocation, (int) SystemComponents.Moleculesize[SelectedComponent]);
210 if(!FF.noCharges && SystemComponents.hasPartialCharge[SelectedComponent])
211 {
212 Update_Ewald_Vector(Sims.Box, false, SystemComponents, SelectedComponent);
213 }
214 Update_NumberOfMolecules(SystemComponents, Sims.d_a, SelectedComponent, INSERTION);
215 return energy;
216 }
217 energy.zero();
218 return energy;
219}
220//Zhao's note: This insertion only takes care of the full (not fractional) molecules//
221static inline MoveEnergy Insertion(Components& SystemComponents, Simulations& Sims, ForceField& FF, RandomNumber& Random, WidomStruct& Widom, size_t SelectedMolInComponent, size_t SelectedComponent)
222{
223 //Get Number of Molecules for this component (For updating TMMC)//
224 //This is the OLD STATE//
225 double NMol = SystemComponents.NumberOfMolecule_for_Component[SelectedComponent];
226 if(SystemComponents.hasfractionalMolecule[SelectedComponent]) NMol--;
227 double TMMCPacc = 0.0;
228
229 SystemComponents.Moves[SelectedComponent].InsertionTotal ++;
230 bool SuccessConstruction = false;
231 double Rosenbluth = 0.0;
232 size_t SelectedTrial = 0;
233 double preFactor = 0.0;
234
235 double2 newScale = SystemComponents.Lambda[SelectedComponent].SET_SCALE(1.0); //Set scale for full molecule (lambda = 1.0)//
236 MoveEnergy energy = Insertion_Body(SystemComponents, Sims, FF, Random, Widom, SelectedMolInComponent, SelectedComponent, Rosenbluth, SuccessConstruction, SelectedTrial, preFactor, false, newScale);
237 if(!SuccessConstruction)
238 {
239 //If unsuccessful move (Overlap), Pacc = 0//
240 SystemComponents.Tmmc[SelectedComponent].Update(TMMCPacc, NMol, INSERTION);
241 SystemComponents.Moves[SelectedComponent].RecordRosen(0.0, INSERTION);
242 energy.zero();
243 return energy;
244 }
245 double IdealRosen = SystemComponents.IdealRosenbluthWeight[SelectedComponent];
246 double RANDOM = Get_Uniform_Random();
247 TMMCPacc = preFactor * Rosenbluth / IdealRosen; //Unbiased Acceptance//
248 //Apply the bias according to the macrostate//
249 SystemComponents.Tmmc[SelectedComponent].ApplyWLBias(preFactor, NMol, INSERTION);
250 SystemComponents.Tmmc[SelectedComponent].ApplyTMBias(preFactor, NMol, INSERTION);
251
252 bool Accept = false;
253 if(RANDOM < preFactor * Rosenbluth / IdealRosen) Accept = true;
254 SystemComponents.Tmmc[SelectedComponent].TreatAccOutofBound(Accept, NMol, INSERTION);
255 SystemComponents.Moves[SelectedComponent].RecordRosen(Rosenbluth, INSERTION);
256
257 if(Accept)
258 { // accept the move
259 SystemComponents.Moves[SelectedComponent].InsertionAccepted ++;
260 AcceptInsertion(SystemComponents, Sims, SelectedComponent, SelectedTrial, FF.noCharges, INSERTION);
261 SystemComponents.Tmmc[SelectedComponent].Update(TMMCPacc, NMol, INSERTION);
262 return energy;
263 }
264 //else
265 //Zhao's note: Even if the move is rejected by acceptance rule, still record the Pacc//
266 SystemComponents.Tmmc[SelectedComponent].Update(TMMCPacc, NMol, INSERTION);
267 energy.zero();
268 return energy;
269}
270
271static inline MoveEnergy Deletion(Components& SystemComponents, Simulations& Sims, ForceField& FF, RandomNumber& Random, WidomStruct& Widom, size_t SelectedMolInComponent, size_t SelectedComponent)
272{
273 //Get Number of Molecules for this component (For updating TMMC)//
274 double NMol = SystemComponents.NumberOfMolecule_for_Component[SelectedComponent];
275 if(SystemComponents.hasfractionalMolecule[SelectedComponent]) NMol--;
276 double TMMCPacc = 0.0;
277
278 SystemComponents.Moves[SelectedComponent].DeletionTotal ++;
279
280 double preFactor = 0.0;
281 bool SuccessConstruction = false;
282 MoveEnergy energy;
283 double Rosenbluth = 0.0;
284 size_t UpdateLocation = 0;
285 int CBMCType = CBMC_DELETION; //Deletion//
286 double2 Scale = SystemComponents.Lambda[SelectedComponent].SET_SCALE(1.0); //Set scale for full molecule (lambda = 1.0), Zhao's note: This is not used in deletion, just set to 1//
287 //Wrapper for the deletion move//
288 energy = Deletion_Body(SystemComponents, Sims, FF, Random, Widom, SelectedMolInComponent, SelectedComponent, UpdateLocation, Rosenbluth, SuccessConstruction, preFactor, Scale);
289 if(!SuccessConstruction)
290 {
291 //If unsuccessful move (Overlap), Pacc = 0//
292 SystemComponents.Tmmc[SelectedComponent].Update(TMMCPacc, NMol, DELETION);
293 SystemComponents.Moves[SelectedComponent].RecordRosen(0.0, DELETION);
294 energy.zero();
295 return energy;
296 }
297 //DEBUG//
298 /*
299 if(SystemComponents.CURRENTCYCLE == 28981)
300 { printf("Selected Molecule: %zu\n", SelectedMolInComponent);
301 printf("DELETION MOVE ENERGY: "); energy.print();
302 }
303 */
304 double IdealRosen = SystemComponents.IdealRosenbluthWeight[SelectedComponent];
305 double RANDOM = Get_Uniform_Random();
306 TMMCPacc = preFactor * IdealRosen / Rosenbluth; //Unbiased Acceptance//
307 //Apply the bias according to the macrostate//
308 SystemComponents.Tmmc[SelectedComponent].ApplyWLBias(preFactor, NMol, DELETION);
309 SystemComponents.Tmmc[SelectedComponent].ApplyTMBias(preFactor, NMol, DELETION);
310
311 bool Accept = false;
312 if(RANDOM < preFactor * IdealRosen / Rosenbluth) Accept = true;
313 SystemComponents.Tmmc[SelectedComponent].TreatAccOutofBound(Accept, NMol, DELETION);
314
315 if(Accept)
316 { // accept the move
317 SystemComponents.Moves[SelectedComponent].DeletionAccepted ++;
318 AcceptDeletion(SystemComponents, Sims, SelectedComponent, UpdateLocation, SelectedMolInComponent, FF.noCharges);
319 //If unsuccessful move (Overlap), Pacc = 0//
320 SystemComponents.Tmmc[SelectedComponent].Update(TMMCPacc, NMol, DELETION);
321 energy.take_negative();
322 return energy;
323 }
324 //Zhao's note: Even if the move is rejected by acceptance rule, still record the Pacc//
325 SystemComponents.Tmmc[SelectedComponent].Update(TMMCPacc, NMol, DELETION);
326 energy.zero();
327 return energy;
328}
329
330static inline void GibbsParticleTransfer(std::vector<Components>& SystemComponents, Simulations*& Sims, ForceField FF, RandomNumber Random, std::vector<WidomStruct>& Widom, std::vector<SystemEnergies>& Energy, size_t SelectedComponent, Gibbs& GibbsStatistics)
331{
332 size_t NBox = SystemComponents.size();
333 size_t SelectedBox = 0;
334 size_t OtherBox = 1;
335 GibbsStatistics.GibbsXferStats.x += 1.0;
336
337 if(Get_Uniform_Random() > 0.5)
338 {
339 SelectedBox = 1;
340 OtherBox = 0;
341 }
342
343 //printf("Performing Gibbs Particle Transfer Move! on Box[%zu]\n", SelectedBox);
344 double2 Scale = {1.0, 1.0};
345 bool TransferFractionalMolecule = false;
346 //Randomly Select a Molecule from the OtherBox in the SelectedComponent//
347 if(SystemComponents[OtherBox].NumberOfMolecule_for_Component[SelectedComponent] == 0) return;
348 size_t DeletionSelectedMol = (size_t) (Get_Uniform_Random()*(SystemComponents[OtherBox].NumberOfMolecule_for_Component[SelectedComponent]));
349 //Special treatment for transfering the fractional molecule//
350 if(SystemComponents[OtherBox].hasfractionalMolecule[SelectedComponent] && SystemComponents[OtherBox].Lambda[SelectedComponent].FractionalMoleculeID == DeletionSelectedMol)
351 {
352 double oldBin = SystemComponents[OtherBox].Lambda[SelectedComponent].currentBin;
353 double delta = SystemComponents[OtherBox].Lambda[SelectedComponent].delta;
354 double oldLambda = delta * static_cast<double>(oldBin);
355 Scale = SystemComponents[OtherBox].Lambda[SelectedComponent].SET_SCALE(oldLambda);
356 TransferFractionalMolecule = true;
357 }
358 //Perform Insertion on the selected System, then deletion on the other system//
359
361 // PERFORMING INSERTION ON THE SELECTED SYSTEM //
363 size_t InsertionSelectedTrial = 0;
364 double InsertionPrefactor = 0.0;
365 double InsertionRosen = 0.0;
366 bool InsertionSuccess = false;
367 size_t InsertionSelectedMol = 0; //It is safer to ALWAYS choose the first atom as the template for CBMC_INSERTION//
368 MoveEnergy InsertionEnergy = Insertion_Body(SystemComponents[SelectedBox], Sims[SelectedBox], FF, Random, Widom[SelectedBox], InsertionSelectedMol, SelectedComponent, InsertionRosen, InsertionSuccess, InsertionSelectedTrial, InsertionPrefactor, false, Scale);
369 printf("Gibbs Particle Insertion energy: "); InsertionEnergy.print();
370 if(!InsertionSuccess) return;
371
373 // PERFORMING DELETION ON THE OTHER SYSTEM //
375 size_t DeletionUpdateLocation = 0;
376 double DeletionPrefactor = 0.0;
377 double DeletionRosen = 0.0;
378 bool DeletionSuccess = false;
379 MoveEnergy DeletionEnergy = Deletion_Body(SystemComponents[OtherBox], Sims[OtherBox], FF, Random, Widom[OtherBox], DeletionSelectedMol, SelectedComponent, DeletionUpdateLocation, DeletionRosen, DeletionSuccess, DeletionPrefactor, Scale);
380 printf("Gibbs Particle Deletion energy: "); DeletionEnergy.print();
381 if(!DeletionSuccess) return;
382
383 bool Accept = false;
384
385 double NMolA= static_cast<double>(SystemComponents[SelectedBox].TotalNumberOfMolecules - SystemComponents[SelectedBox].NumberOfFrameworks);
386 double NMolB= static_cast<double>(SystemComponents[OtherBox].TotalNumberOfMolecules - SystemComponents[OtherBox].NumberOfFrameworks);
387 //Minus the fractional molecule//
388 for(size_t comp = 0; comp < SystemComponents[SelectedBox].Total_Components; comp++)
389 {
390 if(SystemComponents[SelectedBox].hasfractionalMolecule[comp])
391 {
392 NMolA-=1.0;
393 }
394 }
395 for(size_t comp = 0; comp < SystemComponents[OtherBox].Total_Components; comp++)
396 {
397 if(SystemComponents[OtherBox].hasfractionalMolecule[comp])
398 {
399 NMolB-=1.0;
400 }
401 }
402
403 //This assumes that the two boxes share the same temperature, it might not be true//
404 if(Get_Uniform_Random()< (InsertionRosen * NMolB * Sims[SelectedBox].Box.Volume) / (DeletionRosen * (NMolA + 1) * Sims[OtherBox].Box.Volume)) Accept = true;
405
406 //printf("SelectedBox: %zu, OtherBox: %zu, InsertionEnergy: %.5f(%.5f %.5f), DeletionEnergy: %.5f(%.5f %.5f)\n", SelectedBox, OtherBox, InsertionEnergy, SystemComponents[SelectedBox].tempdeltaVDWReal, SystemComponents[SelectedBox].tempdeltaEwald, DeletionEnergy, SystemComponents[OtherBox].tempdeltaVDWReal, SystemComponents[OtherBox].tempdeltaEwald);
407
408 if(Accept)
409 {
410 GibbsStatistics.GibbsXferStats.y += 1.0;
411 // Zhao's note: the assumption for the below implementation is that the component index are the same for both systems //
412 // for example, methane in box A is component 1, it has to be component 1 in box B //
413 //For the box that is deleting the molecule, update the recorded fractional molecule ID//
414 // Update System information regarding the fractional molecule, if the fractional molecule is being transfered //
415 if(TransferFractionalMolecule)
416 {
417 SystemComponents[SelectedBox].hasfractionalMolecule[SelectedComponent] = true;
418 SystemComponents[OtherBox].hasfractionalMolecule[SelectedComponent] = false;
419 SystemComponents[SelectedBox].Lambda[SelectedComponent].currentBin = SystemComponents[OtherBox].Lambda[SelectedComponent].currentBin;
420 }
422 // UPDATE INSERTION FOR THE SELECTED BOX //
424 AcceptInsertion(SystemComponents[SelectedBox], Sims[SelectedBox], SelectedComponent, InsertionSelectedTrial, FF.noCharges, INSERTION);
425
426 SystemComponents[SelectedBox].deltaE += InsertionEnergy;
427 //Energy[SelectedBox].running_energy += InsertionEnergy.total();
428 //printf("Insert Box: %zu, Insertion Energy: %.5f\n", SelectedBox, InsertionEnergy);
429
431 // UPDATE DELETION FOR THE OTHER BOX //
433 AcceptDeletion(SystemComponents[OtherBox], Sims[OtherBox], SelectedComponent, DeletionUpdateLocation, DeletionSelectedMol, FF.noCharges);
434 Energy[OtherBox].running_energy -= DeletionEnergy.total();
435 SystemComponents[OtherBox].deltaE -= DeletionEnergy;
436 //printf("Delete Box: %zu, Insertion Energy: %.5f\n", OtherBox, DeletionEnergy);
437 }
438}
439
440__global__ void copy_firstbead_to_new(Atoms NEW, Atoms* d_a, size_t comp, size_t position)
441{
442 NEW.pos[0] = d_a[comp].pos[position];
443}
444
445__global__ void Update_IdentitySwap_Insertion_data(Atoms* d_a, double3* temp, size_t NEWComponent, size_t UpdateLocation, size_t MolID, size_t Molsize)
446{
447 //Zhao's note: assuming not working with fractional molecule for Identity swap//
448 for(size_t i = 0; i < Molsize; i++)
449 {
450 size_t realLocation = UpdateLocation + i;
451 d_a[NEWComponent].pos[realLocation] = temp[i];
452 d_a[NEWComponent].scale[realLocation] = 1.0;
453 d_a[NEWComponent].charge[realLocation]= d_a[NEWComponent].charge[i];
454 d_a[NEWComponent].scaleCoul[realLocation]=1.0;
455 d_a[NEWComponent].Type[realLocation] = d_a[NEWComponent].Type[i];
456 d_a[NEWComponent].MolID[realLocation] = MolID;
457 }
458 d_a[NEWComponent].size += Molsize;
459}
460
461static inline MoveEnergy IdentitySwapMove(Components& SystemComponents, Simulations& Sims, WidomStruct& Widom, ForceField& FF, RandomNumber& Random)
462{
463 //Identity Swap is sort-of Reinsertion//
464 //The difference is that the component of the molecule are different//
465 //Selected Molecule is the molecule that is being identity-swapped//
466
467 //Zhao's note: If CO2/CH4 mixture, since CH4 doesn't have rotation moves, identity swap will be performed more times for CH4 than for CO2. Add this switch of Old/NewComponent to avoid this issue.
468 size_t NEWComponent = 0;
469 size_t OLDComponent = 0;
470 size_t NOld = 0;
471 //It seems that identity swap can swap back to its own species, so need to relax the current criterion//
472 //while(NEWComponent == OLDComponent || NEWComponent == 0 || NEWComponent >= SystemComponents.NComponents.x)
473 //Must select adsorbate species, cannot exceed number of species in sim, oldcomponent number of mol > 0//
474 if((SystemComponents.TotalNumberOfMolecules - SystemComponents.NumberOfFrameworks) == 0) //No adsorbates
475 {
476 MoveEnergy Empty;
477 return Empty;
478 }
479 while(OLDComponent == 0 || OLDComponent >= SystemComponents.NComponents.x ||
480 NEWComponent == 0 || NEWComponent >= SystemComponents.NComponents.x ||
481 NOld == 0)
482 {
483 OLDComponent = (size_t) (Get_Uniform_Random()*(SystemComponents.NComponents.x - SystemComponents.NComponents.y)) + SystemComponents.NComponents.y;
484 NEWComponent = (size_t) (Get_Uniform_Random()*(SystemComponents.NComponents.x - SystemComponents.NComponents.y)) + SystemComponents.NComponents.y;
485 NOld = SystemComponents.NumberOfMolecule_for_Component[OLDComponent];
486 /*
487 if(Get_Uniform_Random() < 0.5)
488 {
489 size_t tempComponent = OLDComponent;
490 OLDComponent = NEWComponent;
491 NEWComponent = tempComponent;
492 }
493 */
494 }
495
496 size_t OLDMolInComponent = (size_t) (Get_Uniform_Random() * SystemComponents.NumberOfMolecule_for_Component[OLDComponent]);
497
498 //printf("Cycle: %zu, OLDComp: %zu, NEWComp: %zu\n", SystemComponents.CURRENTCYCLE, OLDComponent, NEWComponent);
499 //JUST FOR DEBUG//
500 //NEWComponent = 1; OLDComponent = 1;
501 //
502
503 double NNEWMol = static_cast<double>(SystemComponents.NumberOfMolecule_for_Component[NEWComponent]);
504 double NOLDMol = static_cast<double>(SystemComponents.NumberOfMolecule_for_Component[OLDComponent]);
505 if(SystemComponents.hasfractionalMolecule[NEWComponent]) NNEWMol -= 1.0;
506 if(SystemComponents.hasfractionalMolecule[OLDComponent]) NOLDMol -= 1.0;
507
508 //FOR DEBUG//
509 //OLDMolInComponent = 52;
510
511 size_t NEWMolInComponent = SystemComponents.NumberOfMolecule_for_Component[NEWComponent];
512
513 SystemComponents.Moves[OLDComponent].IdentitySwapRemoveTotal ++;
514 SystemComponents.Moves[NEWComponent].IdentitySwapAddTotal ++;
515 SystemComponents.Moves[OLDComponent].IdentitySwap_Total_TO[NEWComponent] ++;
516
517 bool SuccessConstruction = false;
518 size_t SelectedTrial = 0;
519 MoveEnergy energy; MoveEnergy old_energy; double StoredR = 0.0;
520
522 // INSERTION //
524 int CBMCType = IDENTITY_SWAP_NEW; //Reinsertion-Insertion//
525 //Take care of the frist bead location HERE//
526 copy_firstbead_to_new<<<1,1>>>(Sims.New, Sims.d_a, OLDComponent, OLDMolInComponent * SystemComponents.Moleculesize[OLDComponent]);
527 //WRITE THE component and molecule ID THAT IS GOING TO BE IGNORED DURING THE ENERGY CALC//
528 Sims.ExcludeList[0] = {OLDComponent, OLDMolInComponent};
529 double2 newScale = SystemComponents.Lambda[NEWComponent].SET_SCALE(1.0); //Zhao's note: not used in reinsertion, just set to 1.0//
530 double Rosenbluth=Widom_Move_FirstBead_PARTIAL(SystemComponents, Sims, FF, Random, Widom, NEWMolInComponent, NEWComponent, CBMCType, StoredR, &SelectedTrial, &SuccessConstruction, &energy, newScale);
531
532 if(Rosenbluth <= 1e-150) SuccessConstruction = false; //Zhao's note: added this protection bc of weird error when testing GibbsParticleXfer
533
534 if(!SuccessConstruction)
535 {
536 energy.zero(); Sims.ExcludeList[0] = {-1, -1}; //Set to negative so that excludelist is ignored
537 return energy;
538 }
539
540 if(SystemComponents.Moleculesize[NEWComponent] > 1 && Rosenbluth > 1e-150)
541 {
542 size_t SelectedFirstBeadTrial = SelectedTrial;
543 MoveEnergy temp_energy = energy;
544 Rosenbluth*=Widom_Move_Chain_PARTIAL(SystemComponents, Sims, FF, Random, Widom, NEWMolInComponent, NEWComponent, CBMCType, &SelectedTrial, &SuccessConstruction, &energy, SelectedFirstBeadTrial, newScale); //True for doing insertion for reinsertion, different in MoleculeID//
545 if(Rosenbluth <= 1e-150) SuccessConstruction = false; //Zhao's note: added this protection bc of weird error when testing GibbsParticleXfer
546 if(!SuccessConstruction)
547 {
548 energy.zero(); Sims.ExcludeList[0] = {-1, -1}; //Set to negative so that excludelist is ignored
549 return energy;
550 }
551 energy += temp_energy;
552 }
553 //printf("NEW MOLECULE ENERGY:"); energy.print();
554
555 // Store The New Locations //
556 double3* temp;
557 cudaMalloc(&temp, sizeof(double3) * SystemComponents.Moleculesize[NEWComponent]);
558 StoreNewLocation_Reinsertion<<<1,1>>>(Sims.Old, Sims.New, temp, SelectedTrial, SystemComponents.Moleculesize[NEWComponent]);
559
561 // RETRACE //
563
564 Sims.ExcludeList[0] = {-1, -1}; //Set to negative so that excludelist is ignored
565
566 CBMCType = IDENTITY_SWAP_OLD; //Identity_Swap for the old configuration//
567 double Old_Rosen=Widom_Move_FirstBead_PARTIAL(SystemComponents, Sims, FF, Random, Widom, OLDMolInComponent, OLDComponent, CBMCType, StoredR, &SelectedTrial, &SuccessConstruction, &old_energy, newScale);
568 if(SystemComponents.Moleculesize[OLDComponent] > 1)
569 {
570 size_t SelectedFirstBeadTrial = SelectedTrial;
571 MoveEnergy temp_energy = old_energy;
572 Old_Rosen*=Widom_Move_Chain_PARTIAL(SystemComponents, Sims, FF, Random, Widom, OLDMolInComponent, OLDComponent, CBMCType, &SelectedTrial, &SuccessConstruction, &old_energy, SelectedFirstBeadTrial, newScale);
573 old_energy += temp_energy;
574 }
575
576 //printf("OLD MOLECULE ENERGY:"); old_energy.print();
577 energy -= old_energy;
578
579 //Calculate Ewald//
580 size_t UpdateLocation = SystemComponents.Moleculesize[OLDComponent] * OLDMolInComponent;
581 if(!FF.noCharges)
582 {
583 double2 EwaldE = GPU_EwaldDifference_IdentitySwap(Sims.Box, Sims.d_a, Sims.Old, temp, FF, Sims.Blocksum, SystemComponents, OLDComponent, NEWComponent, UpdateLocation);
584 energy.GGEwaldE = EwaldE.x;
585 energy.HGEwaldE = EwaldE.y;
586 Rosenbluth *= std::exp(-SystemComponents.Beta * (EwaldE.x + EwaldE.y));
587 //printf("Ewald PreFactor: %.5f\n", ::exp(-SystemComponents.Beta * (EwaldE.x + EwaldE.y)));
588 }
589
590 energy.TailE = TailCorrectionIdentitySwap(SystemComponents, NEWComponent, OLDComponent, FF.size, Sims.Box.Volume);
591 Rosenbluth *= std::exp(-SystemComponents.Beta * energy.TailE);
592
593 //NO DEEP POTENTIAL FOR THIS MOVE!//
594 if(SystemComponents.UseDNNforHostGuest) throw std::runtime_error("NO DEEP POTENTIAL FOR IDENTITY SWAP!");
595
597 // PREPARE ACCEPTANCE RULE //
599 double NEWIdealRosen = SystemComponents.IdealRosenbluthWeight[NEWComponent];
600 double OLDIdealRosen = SystemComponents.IdealRosenbluthWeight[OLDComponent];
601
602 double preFactor = GetPrefactor(SystemComponents, Sims, NEWComponent, INSERTION);
603 preFactor *= GetPrefactor(SystemComponents, Sims, OLDComponent, DELETION);
604
605 double Pacc = preFactor * (Rosenbluth / NEWIdealRosen) / (Old_Rosen / OLDIdealRosen);
606
607 //printf("Rosenbluth Weights %.5f (New), %.5f (Old)\n", Rosenbluth, Old_Rosen);
608 //printf("Partial Fugacities: %.5f (New), %.5f (Old)\n", PartialFugacityNew, PartialFugacityOld);
609 //printf("NNEW: %.5f, NOLD: %.5f\n", NNEWMol, NOLDMol);
610 //printf("PAcc: %.5f\n", Pacc);
611
612 //Determine whether to accept or reject the insertion
613 double RANDOM = Get_Uniform_Random();
614 bool Accept = false;
615 if(RANDOM < Pacc) Accept = true;
616
617 if(Accept)
618 { // accept the move
619 SystemComponents.Moves[NEWComponent].IdentitySwapAddAccepted ++;
620 SystemComponents.Moves[OLDComponent].IdentitySwapRemoveAccepted ++;
621
622 SystemComponents.Moves[OLDComponent].IdentitySwap_Acc_TO[NEWComponent] ++;
623 if(NEWComponent != OLDComponent)
624 {
625 size_t LastMolecule = SystemComponents.NumberOfMolecule_for_Component[OLDComponent]-1;
626 size_t LastLocation = LastMolecule*SystemComponents.Moleculesize[OLDComponent];
627 Update_deletion_data<<<1,1>>>(Sims.d_a, OLDComponent, UpdateLocation, (int) SystemComponents.Moleculesize[OLDComponent], LastLocation);
628
629 //The function below will only be processed if the system has a fractional molecule and the transfered molecule is NOT the fractional one //
630 if((SystemComponents.hasfractionalMolecule[OLDComponent])&&(LastMolecule == SystemComponents.Lambda[OLDComponent].FractionalMoleculeID))
631 {
632 //Since the fractional molecule is moved to the place of the selected deleted molecule, update fractional molecule ID on host
633 SystemComponents.Lambda[OLDComponent].FractionalMoleculeID = OLDMolInComponent;
634 }
635
636 UpdateLocation = SystemComponents.Moleculesize[NEWComponent] * NEWMolInComponent;
637 Update_IdentitySwap_Insertion_data<<<1,1>>>(Sims.d_a, temp, NEWComponent, UpdateLocation, NEWMolInComponent, SystemComponents.Moleculesize[NEWComponent]); checkCUDAError("error Updating Identity Swap Insertion data");
638
639 Update_NumberOfMolecules(SystemComponents, Sims.d_a, NEWComponent, INSERTION);
640 Update_NumberOfMolecules(SystemComponents, Sims.d_a, OLDComponent, DELETION);
641 }
642 else //If they are the same species, just update in the reinsertion way (just the locations)//
643 {
644 //Regarding the UpdateLocation, in the code it is the old molecule position//
645 //if the OLDComponent = NEWComponent, The new location will be filled into the old position//
646 //So just use what is already in the code//
647 Update_Reinsertion_data<<<1,SystemComponents.Moleculesize[OLDComponent]>>>(Sims.d_a, temp, OLDComponent, UpdateLocation);
648 }
649 cudaFree(temp);
650 //Zhao's note: BUG!!!!, Think about if OLD/NEW Component belong to different type (framework/adsorbate)//
651 if(!FF.noCharges && ((SystemComponents.hasPartialCharge[NEWComponent]) ||(SystemComponents.hasPartialCharge[OLDComponent])))
652 Update_Ewald_Vector(Sims.Box, false, SystemComponents, NEWComponent);
653 //energy.print();
654 return energy;
655 }
656 else
657 {
658 cudaFree(temp);
659 energy.zero();
660 return energy;
661 }
662}
Definition data_struct.h:746
Definition data_struct.h:843
Definition data_struct.h:794
Definition data_struct.h:61
Definition data_struct.h:615
Definition data_struct.h:1053
Definition data_struct.h:1027
Definition data_struct.h:1044