My Project
Loading...
Searching...
No Matches
mc_swap_utilities.h
1static inline MoveEnergy Insertion_Body(Components& SystemComponents, Simulations& Sims, ForceField& FF, RandomNumber& Random, WidomStruct& Widom, size_t SelectedMolInComponent, size_t SelectedComponent, double& Rosenbluth, bool& SuccessConstruction, size_t& SelectedTrial, double& preFactor, bool previous_step, double2 newScale)
2{
3 MoveEnergy energy; double StoredR = 0.0;
4 int CBMCType = CBMC_INSERTION; //Insertion//
5 Rosenbluth=Widom_Move_FirstBead_PARTIAL(SystemComponents, Sims, FF, Random, Widom, SelectedMolInComponent, SelectedComponent, CBMCType, StoredR, &SelectedTrial, &SuccessConstruction, &energy, newScale); //Not reinsertion, not Retrace//
6
7 if(Rosenbluth <= 1e-150) SuccessConstruction = false; //Zhao's note: added this protection bc of weird error when testing GibbsParticleXfer
8 if(!SuccessConstruction)
9 {
10 //printf("Early return FirstBead\n");
11 energy.zero();
12 return energy;
13 }
14 if(SystemComponents.Moleculesize[SelectedComponent] > 1 && Rosenbluth > 1e-150)
15 {
16 size_t SelectedFirstBeadTrial = SelectedTrial; MoveEnergy temp_energy = energy;
17 Rosenbluth*=Widom_Move_Chain_PARTIAL(SystemComponents, Sims, FF, Random, Widom, SelectedMolInComponent, SelectedComponent, CBMCType, &SelectedTrial, &SuccessConstruction, &energy, SelectedFirstBeadTrial, newScale);
18
19 if(Rosenbluth <= 1e-150) SuccessConstruction = false;
20 if(!SuccessConstruction)
21 {
22 //printf("Early return Chain\n");
23 energy.zero();
24 return energy;
25 }
26 energy += temp_energy;
27 }
28 //Determine whether to accept or reject the insertion
29 preFactor = GetPrefactor(SystemComponents, Sims, SelectedComponent, INSERTION);
30 //Ewald Correction, done on HOST (CPU) //
31 int MoveType = INSERTION; //Normal Insertion, including fractional insertion, no previous step (do not use temprorary tempEik)//
32 if(previous_step) //Fractional Insertion after a lambda change move that makes the old fractional molecule full//
33 {
34 MoveType = CBCF_INSERTION; // CBCF fractional insertion //
35 }
36 bool EwaldPerformed = false;
37 if(!FF.noCharges && SystemComponents.hasPartialCharge[SelectedComponent])
38 {
39 double2 EwaldE = GPU_EwaldDifference_General(Sims.Box, Sims.d_a, Sims.New, Sims.Old, FF, Sims.Blocksum, SystemComponents, SelectedComponent, MoveType, SelectedTrial, newScale);
40
41 energy.GGEwaldE = EwaldE.x;
42 energy.HGEwaldE = EwaldE.y;
43 Rosenbluth *= std::exp(-SystemComponents.Beta * (EwaldE.x + EwaldE.y));
44 EwaldPerformed = true;
45 }
46 double TailE = 0.0;
47 TailE = TailCorrectionDifference(SystemComponents, SelectedComponent, FF.size, Sims.Box.Volume, MoveType);
48 Rosenbluth *= std::exp(-SystemComponents.Beta * TailE);
49 energy.TailE= TailE;
50 //Calculate DNN//
51 //Put it after Ewald summation, the required positions are already in place (done by the preparation parts of Ewald Summation)//
52 if(SystemComponents.UseDNNforHostGuest)
53 {
54 if(!EwaldPerformed) Prepare_DNN_InitialPositions(Sims.d_a, Sims.New, Sims.Old, SystemComponents, SelectedComponent, MoveType, SelectedTrial);
55 double DNN_New = DNN_Prediction_Move(SystemComponents, Sims, SelectedComponent, INSERTION);
56 energy.DNN_E = DNN_New;
57 double correction = energy.DNN_Correction();
58 if(fabs(correction) > SystemComponents.DNNDrift) //If there is a huge drift in the energy correction between DNN and Classical HostGuest//
59 {
60 //printf("INSERTION: Bad Prediction, reject the move!!!\n");
61 SystemComponents.InsertionDNNReject ++;
62 SuccessConstruction = false;
63 WriteOutliers(SystemComponents, Sims, DNN_INSERTION, energy, correction);
64 energy.zero();
65 return energy;
66 }
67 SystemComponents.InsertionDNNDrift += fabs(correction);
68 Rosenbluth *= std::exp(-SystemComponents.Beta * correction);
69 }
70 //printf("Insertion energy summary: "); energy.print();
71 return energy;
72}
73
74static inline MoveEnergy Deletion_Body(Components& SystemComponents, Simulations& Sims, ForceField& FF, RandomNumber& Random, WidomStruct& Widom, size_t SelectedMolInComponent, size_t SelectedComponent, size_t& UpdateLocation, double& Rosenbluth, bool& SuccessConstruction, double& preFactor, double2 Scale)
75{
76 size_t SelectedTrial = 0;
77 MoveEnergy energy;
78
79 double StoredR = 0.0; //Don't use this for Deletion//
80 int CBMCType = CBMC_DELETION; //Deletion//
81 //Zhao's note: Deletion_body will be part the GibbsParticleTransferMove, and the fractional molecule might be selected, so Scale will not be 1.0//
82 //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//
83 Rosenbluth=Widom_Move_FirstBead_PARTIAL(SystemComponents, Sims, FF, Random, Widom, SelectedMolInComponent, SelectedComponent, CBMCType, StoredR, &SelectedTrial, &SuccessConstruction, &energy, Scale);
84 if(Rosenbluth <= 1e-150) SuccessConstruction = false;
85 if(!SuccessConstruction)
86 {
87 energy.zero();
88 return energy;
89 }
90
91 //DEBUG//
92 /*
93 if(SystemComponents.CURRENTCYCLE == 28981)
94 {
95 printf("DELETION FIRST BEAD ENERGY: "); energy.print();
96 printf("EXCLUSION LIST: %d %d\n", Sims.ExcludeList[0].x, Sims.ExcludeList[0].y);
97 }
98 */
99
100 if(SystemComponents.Moleculesize[SelectedComponent] > 1)
101 {
102 size_t SelectedFirstBeadTrial = SelectedTrial; MoveEnergy temp_energy = energy;
103 Rosenbluth*=Widom_Move_Chain_PARTIAL(SystemComponents, Sims, FF, Random, Widom, SelectedMolInComponent, SelectedComponent, CBMCType, &SelectedTrial, &SuccessConstruction, &energy, SelectedFirstBeadTrial, Scale); //The false is for Reinsertion//
104 energy += temp_energy;
105 }
106 if(Rosenbluth <= 1e-150) SuccessConstruction = false;
107 if(!SuccessConstruction)
108 {
109 energy.zero();
110 return energy;
111 }
112
113 preFactor = GetPrefactor(SystemComponents, Sims, SelectedComponent, DELETION);
114
115 UpdateLocation = SelectedMolInComponent * SystemComponents.Moleculesize[SelectedComponent];
116 bool EwaldPerformed = false;
117 if(!FF.noCharges && SystemComponents.hasPartialCharge[SelectedComponent])
118 {
119 int MoveType = DELETION; if(Scale.y < 1.0) MoveType = CBCF_DELETION;
120
121 double2 EwaldE = GPU_EwaldDifference_General(Sims.Box, Sims.d_a, Sims.New, Sims.Old, FF, Sims.Blocksum, SystemComponents, SelectedComponent, MoveType, UpdateLocation, Scale);
122 Rosenbluth /= std::exp(-SystemComponents.Beta * (EwaldE.x + EwaldE.y));
123 energy.GGEwaldE = -1.0 * EwaldE.x;
124 energy.HGEwaldE = -1.0 * EwaldE.y; //Becareful with the sign here, you need a HG sum, but HGVDWReal and HGEwaldE here have opposite signs???//
125 EwaldPerformed = true;
126 }
127 double TailE = 0.0;
128 TailE = TailCorrectionDifference(SystemComponents, SelectedComponent, FF.size, Sims.Box.Volume, DELETION);
129 Rosenbluth /= std::exp(-SystemComponents.Beta * TailE);
130 energy.TailE = -TailE;
131 //Calculate DNN//
132 //Put it after Ewald summation, the required positions are already in place (done by the preparation parts of Ewald Summation)//
133 if(SystemComponents.UseDNNforHostGuest)
134 {
135 if(!EwaldPerformed) Prepare_DNN_InitialPositions(Sims.d_a, Sims.New, Sims.Old, SystemComponents, SelectedComponent, DELETION, UpdateLocation);
136 //Deletion positions stored in Old//
137 double DNN_New = DNN_Prediction_Move(SystemComponents, Sims, SelectedComponent, DELETION);
138 energy.DNN_E = DNN_New;
139 double correction = energy.DNN_Correction();
140 if(fabs(correction) > SystemComponents.DNNDrift) //If there is a huge drift in the energy correction between DNN and Classical HostGuest//
141 {
142 //printf("DELETION: Bad Prediction, reject the move!!!\n");
143 SystemComponents.DeletionDNNReject ++;
144 SuccessConstruction = false;
145 WriteOutliers(SystemComponents, Sims, DNN_DELETION, energy, correction);
146 energy.zero();
147 return energy;
148 }
149 SystemComponents.DeletionDNNDrift += fabs(correction);
150 Rosenbluth *= std::exp(-SystemComponents.Beta * correction);
151 }
152 return energy;
153}
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