My Project
Loading...
Searching...
No Matches
mc_single_particle.h
1#include "mc_utilities.h"
2
4// Generalized function for single Body moves //
6static inline MoveEnergy SingleBodyMove(Components& SystemComponents, Simulations& Sims, WidomStruct& Widom, ForceField& FF, RandomNumber& Random, size_t SelectedMolInComponent, size_t SelectedComponent, int MoveType)
7{
8 //Get Number of Molecules for this component (For updating TMMC)//
9 double NMol = SystemComponents.NumberOfMolecule_for_Component[SelectedComponent];
10 if(SystemComponents.hasfractionalMolecule[SelectedComponent]) NMol--;
11
12 bool Do_New = false;
13 bool Do_Old = false;
14
15 size_t Atomsize = 0;
16 for(size_t ijk = 0; ijk < SystemComponents.Total_Components; ijk++)
17 Atomsize += SystemComponents.Moleculesize[ijk] * SystemComponents.NumberOfMolecule_for_Component[ijk];
18 size_t Molsize = SystemComponents.Moleculesize[SelectedComponent]; //Get the size of the selected Molecule
19 //Set up Old position and New position arrays
20 if(Molsize >= 1024)
21 {
22 throw std::runtime_error("Molecule size is greater than allocated size, Why so big?\n");
23 }
24 size_t start_position = SelectedMolInComponent*SystemComponents.Moleculesize[SelectedComponent];
25
26 SystemComponents.Moves[SelectedComponent].Record_Move_Total(MoveType);
27 double3 MaxChange = {0.0, 0.0, 0.0};
28 bool CheckOverlap = true;
29 switch (MoveType)
30 {
31 case TRANSLATION:
32 {
33 Do_New = true; Do_Old = true;
34 MaxChange = SystemComponents.MaxTranslation[SelectedComponent];
35 break;
36 }
37 case ROTATION:
38 {
39 Do_New = true; Do_Old = true;
40 MaxChange = SystemComponents.MaxRotation[SelectedComponent];
41 break;
42 }
43 case SPECIAL_ROTATION:
44 {
45 Do_New = true; Do_Old = true;
46 MaxChange = SystemComponents.MaxSpecialRotation[SelectedComponent];
47 //Zhao's note: if we separate framework components, there might be lots of overlaps between different species (node and linker overlaps), we can turn this Overlap flag off//
48 CheckOverlap = false;
49 //printf("Performing move on %zu comp, %zu mol\n", SelectedComponent, SelectedMolInComponent);
50 break;
51 }
52 case SINGLE_INSERTION:
53 {
54 Do_New = true;
55 start_position = 0;
56 break;
57 }
58 case SINGLE_DELETION:
59 {
60 Do_Old = true;
61 break;
62 }
63 }
64 if(!Do_New && !Do_Old) throw std::runtime_error("Doing Nothing For Single Particle Move?\n");
65
66 //Zhao's note: possible bug, you may only need 3 instead of 3 * N random numbers//
67 Random.Check(Molsize);
68 get_new_position<<<1, Molsize>>>(Sims, FF, start_position, SelectedComponent, MaxChange, Random.device_random, Random.offset, MoveType);
69 Random.Update(Molsize);
70
71 // Setup for the pairwise calculation //
72 // New Features: divide the Blocks into two parts: Host-Guest + Guest-Guest //
73
74 size_t NHostAtom = 0; size_t NGuestAtom = 0;
75 for(size_t i = 0; i < SystemComponents.NComponents.y; i++)
76 NHostAtom += SystemComponents.Moleculesize[i] * SystemComponents.NumberOfMolecule_for_Component[i];
77 for(size_t i = SystemComponents.NComponents.y; i < SystemComponents.NComponents.x; i++)
78 NGuestAtom+= SystemComponents.Moleculesize[i] * SystemComponents.NumberOfMolecule_for_Component[i];
79
80 //Zhao's note: Cross term, if the selected species is host atom, the crossAtom = guest, vice versa//
81 size_t NCrossAtom = NHostAtom;
82 if(SelectedComponent < SystemComponents.NComponents.y) //Framework component//
83 NCrossAtom = NGuestAtom;
84 size_t HH_Nthread=0; size_t HH_Nblock=0; Setup_threadblock(NHostAtom * Molsize, &HH_Nblock, &HH_Nthread);
85 size_t HG_Nthread=0; size_t HG_Nblock=0; Setup_threadblock(NCrossAtom * Molsize, &HG_Nblock, &HG_Nthread);
86 size_t GG_Nthread=0; size_t GG_Nblock=0; Setup_threadblock(NGuestAtom * Molsize, &GG_Nblock, &GG_Nthread);
87
88 size_t SameTypeNthread = 0;
89 if(SelectedComponent < SystemComponents.NComponents.y) //Framework-Framework + Framework-Adsorbate//
90 {GG_Nthread = 0; GG_Nblock = 0; SameTypeNthread = HH_Nthread; }
91 else //Framework-Adsorbate + Adsorbate-Adsorbate//
92 {HH_Nthread = 0; HH_Nblock = 0; SameTypeNthread = GG_Nthread; }
93
94 size_t Nthread = std::max(SameTypeNthread, HG_Nthread);
95 size_t Total_Nblock = HH_Nblock + HG_Nblock + GG_Nblock;
96
97 int3 NBlocks = {(int) HH_Nblock, (int) HG_Nblock, (int) GG_Nblock}; //x: HH_Nblock, y: HG_Nblock, z: GG_Nblock;
98 //printf("Total_Comp: %zu, Host Comp: %zu, Adsorbate Comp: %zu\n", SystemComponents.NComponents.x, SystemComponents.NComponents.y, SystemComponents.NComponents.z);
99 //printf("NHostAtom: %zu, HH_Nblock: %zu, HG_Nblock: %zu, NGuestAtom: %zu, GG_Nblock: %zu\n", NHostAtom, HH_Nblock, HG_Nblock, NGuestAtom, GG_Nblock);
100 if(Atomsize != 0)
101 {
102 Calculate_Single_Body_Energy_SEPARATE_HostGuest_VDWReal<<<Total_Nblock, Nthread, Nthread * 2 * sizeof(double)>>>(Sims.Box, Sims.d_a, Sims.Old, Sims.New, FF, Sims.Blocksum, SelectedComponent, Atomsize, Molsize, Sims.device_flag, NBlocks, Do_New, Do_Old, SystemComponents.NComponents);
103
104 cudaMemcpy(SystemComponents.flag, Sims.device_flag, sizeof(bool), cudaMemcpyDeviceToHost);
105 }
106 MoveEnergy tot; bool Accept = false; double Pacc = 0.0;
107 if(!SystemComponents.flag[0] || !CheckOverlap)
108 {
109 double BlockResult[Total_Nblock + Total_Nblock];
110 cudaMemcpy(BlockResult, Sims.Blocksum, 2 * Total_Nblock * sizeof(double), cudaMemcpyDeviceToHost);
111
112 //VDW Part and Real Part Coulomb//
113 for(size_t i = 0; i < HH_Nblock; i++)
114 {
115 tot.HHVDW += BlockResult[i];
116 tot.HHReal+= BlockResult[i + Total_Nblock];
117 //if(MoveType == SPECIAL_ROTATION) printf("HH Block %zu, VDW: %.5f, Real: %.5f\n", i, BlockResult[i], BlockResult[i + Total_Nblock]);
118 }
119 for(size_t i = HH_Nblock; i < HH_Nblock + HG_Nblock; i++)
120 {
121 tot.HGVDW += BlockResult[i];
122 tot.HGReal+= BlockResult[i + Total_Nblock];
123 //printf("HG Block %zu, VDW: %.5f, Real: %.5f\n", i, BlockResult[i], BlockResult[i + Total_Nblock]);
124 }
125 for(size_t i = HH_Nblock + HG_Nblock; i < Total_Nblock; i++)
126 {
127 tot.GGVDW += BlockResult[i];
128 tot.GGReal+= BlockResult[i + Total_Nblock];
129 //printf("GG Block %zu, VDW: %.5f, Real: %.5f\n", i, BlockResult[i], BlockResult[i + Total_Nblock]);
130 }
131
132 /*
133 printf("HG_NBlock: %zu\n", Total_Nblock);
134 printf("Separated VDW : %.5f (HH), %.5f (HG), %.5f (GG)\n", tot.HHVDW, tot.HGVDW , tot.GGVDW);
135 printf("Separated Real: %.5f (HH), %.5f (HG), %.5f (GG)\n", tot.HHReal, tot.HGReal, tot.GGReal);
136 */
137
138 // Calculate Ewald //
139 bool EwaldPerformed = false;
140 if(!FF.noCharges && SystemComponents.hasPartialCharge[SelectedComponent])
141 {
142 double2 newScale = SystemComponents.Lambda[SelectedComponent].SET_SCALE(1.0);
143 double2 EwaldE = GPU_EwaldDifference_General(Sims.Box, Sims.d_a, Sims.New, Sims.Old, FF, Sims.Blocksum, SystemComponents, SelectedComponent, MoveType, 0, newScale);
144 if(HH_Nblock == 0)
145 {
146 tot.GGEwaldE = EwaldE.x;
147 tot.HGEwaldE = EwaldE.y;
148 }
149 else
150 {
151 tot.HHEwaldE = EwaldE.x;
152 tot.HGEwaldE = EwaldE.y;
153 //printf("HHEwald: %.5f, HGEwald: %.5f\n", tot.HHEwaldE, tot.HGEwaldE);
154 }
155 EwaldPerformed = true;
156 }
157 if(SystemComponents.UseDNNforHostGuest)
158 {
159 //Calculate DNN//
160 if(!EwaldPerformed) Prepare_DNN_InitialPositions(Sims.d_a, Sims.New, Sims.Old, SystemComponents, SelectedComponent, MoveType, 0);
161 tot.DNN_E = DNN_Prediction_Move(SystemComponents, Sims, SelectedComponent, MoveType);
162 double correction = tot.DNN_Correction(); //If use DNN, HGVDWReal and HGEwaldE are zeroed//
163 if(fabs(correction) > SystemComponents.DNNDrift) //If there is a huge drift in the energy correction between DNN and Classical HostGuest//
164 {
165 //printf("TRANSLATION/ROTATION: Bad Prediction, reject the move!!!\n");
166 switch(MoveType)
167 {
168 case TRANSLATION: case ROTATION:
169 {
170 SystemComponents.TranslationRotationDNNReject ++; break;
171 }
172 case SINGLE_INSERTION: case SINGLE_DELETION:
173 {
174 SystemComponents.SingleSwapDNNReject ++; break;
175 }
176 }
177 WriteOutliers(SystemComponents, Sims, NEW, tot, correction); //Print New Locations//
178 WriteOutliers(SystemComponents, Sims, OLD, tot, correction); //Print Old Locations//
179 tot.zero();
180 return tot;
181 }
182 SystemComponents.SingleMoveDNNDrift += fabs(correction);
183 }
184
185 double preFactor = GetPrefactor(SystemComponents, Sims, SelectedComponent, MoveType);
186 Pacc = preFactor * std::exp(-SystemComponents.Beta * tot.total());
187
188 //Apply the bias according to the macrostate//
189 if(MoveType == SINGLE_INSERTION || MoveType == SINGLE_DELETION)
190 {
191 SystemComponents.Tmmc[SelectedComponent].ApplyWLBias(preFactor, NMol, MoveType);
192 SystemComponents.Tmmc[SelectedComponent].ApplyTMBias(preFactor, NMol, MoveType);
193 }
194 //if(MoveType == SINGLE_INSERTION) printf("SINGLE INSERTION, tot: %.5f, preFactor: %.5f, Pacc: %.5f\n", tot.total(), preFactor, Pacc);
195 //if(MoveType == SINGLE_DELETION) printf("SINGLE DELETION, tot: %.5f, preFactor: %.5f, Pacc: %.5f\n", tot.total(), preFactor, Pacc);
196 if(Get_Uniform_Random() < preFactor * std::exp(-SystemComponents.Beta * tot.total())) Accept = true;
197 }
198
199 //if(MoveType == SPECIAL_ROTATION)
200 // printf("Framework Component Move, %zu cycle, Molecule: %zu, Energy: %.5f, %s\n", SystemComponents.CURRENTCYCLE, SelectedMolInComponent, tot.total(), Accept ? "Accept" : "Reject");
201 //if(MoveType == SPECIAL_ROTATION) Accept = true;
202
203 switch(MoveType)
204 {
205 case TRANSLATION: case ROTATION: case SPECIAL_ROTATION:
206 {
207 if(Accept)
208 {
209 update_translation_position<<<1,Molsize>>>(Sims.d_a, Sims.New, start_position, SelectedComponent);
210 SystemComponents.Moves[SelectedComponent].Record_Move_Accept(MoveType);
211 if(!FF.noCharges && SystemComponents.hasPartialCharge[SelectedComponent])
212 {
213 Update_Ewald_Vector(Sims.Box, false, SystemComponents, SelectedComponent);
214 }
215 }
216 else {tot.zero(); };
217 SystemComponents.Tmmc[SelectedComponent].Update(1.0, NMol, MoveType);
218 break;
219 }
220 case SINGLE_INSERTION:
221 {
222 SystemComponents.Tmmc[SelectedComponent].TreatAccOutofBound(Accept, NMol, MoveType);
223 if(Accept)
224 {
225 SystemComponents.Moves[SelectedComponent].Record_Move_Accept(MoveType);
226 AcceptInsertion(SystemComponents, Sims, SelectedComponent, 0, FF.noCharges, SINGLE_INSERTION); //0: selectedTrial//
227 }
228 else {tot.zero(); };
229 SystemComponents.Tmmc[SelectedComponent].Update(Pacc, NMol, MoveType);
230 break;
231 }
232 case SINGLE_DELETION:
233 {
234 SystemComponents.Tmmc[SelectedComponent].TreatAccOutofBound(Accept, NMol, MoveType);
235 if(Accept)
236 {
237 SystemComponents.Moves[SelectedComponent].Record_Move_Accept(MoveType);
238 size_t UpdateLocation = SelectedMolInComponent * SystemComponents.Moleculesize[SelectedComponent];
239 AcceptDeletion(SystemComponents, Sims, SelectedComponent, UpdateLocation, SelectedMolInComponent, FF.noCharges);
240 }
241 else {tot.zero(); };
242 SystemComponents.Tmmc[SelectedComponent].Update(Pacc, NMol, MoveType);
243 break;
244 }
245 }
246 //if(MoveType == SINGLE_INSERTION) {printf("Cycle %zu, ENERGY: ", SystemComponents.CURRENTCYCLE); tot.print();}
247 return tot;
248}
Definition data_struct.h:843
Definition data_struct.h:794
Definition data_struct.h:615
Definition data_struct.h:1053
Definition data_struct.h:1027
Definition data_struct.h:1044