1 |
|
#include "Integrator.hpp" |
2 |
|
#include "simError.h" |
3 |
< |
|
3 |
> |
#include <cmath> |
4 |
|
template<typename T> ZConstraint<T>::ZConstraint(SimInfo* theInfo, ForceFields* the_ff) |
5 |
< |
: T(theInfo, the_ff), fz(NULL), indexOfZConsMols(NULL) |
5 |
> |
: T(theInfo, the_ff), fz(NULL), curZPos(NULL), fzOut(NULL), |
6 |
> |
indexOfZConsMols(NULL), forcePolicy(NULL), curZconsTime(0) |
7 |
|
{ |
8 |
|
|
9 |
|
//get properties from SimInfo |
10 |
|
GenericData* data; |
11 |
< |
IndexData* index; |
11 |
> |
ZConsParaData* zConsParaData; |
12 |
|
DoubleData* sampleTime; |
13 |
+ |
DoubleData* tolerance; |
14 |
+ |
StringData* policy; |
15 |
|
StringData* filename; |
16 |
+ |
double COM[3]; |
17 |
+ |
|
18 |
+ |
//by default, the direction of constraint is z |
19 |
+ |
// 0 --> x |
20 |
+ |
// 1 --> y |
21 |
+ |
// 2 --> z |
22 |
+ |
whichDirection = 2; |
23 |
+ |
|
24 |
+ |
//estimate the force constant of harmonical potential |
25 |
+ |
double Kb = 1.986E-3 ; //in kcal/K |
26 |
|
|
27 |
< |
|
28 |
< |
data = info->getProperty("zconsindex"); |
16 |
< |
if(!data) { |
27 |
> |
double halfOfLargestBox = max(info->boxL[0], max(info->boxL[1], info->boxL[2])) /2; |
28 |
> |
zForceConst = Kb * info->target_temp /(halfOfLargestBox * halfOfLargestBox); |
29 |
|
|
30 |
+ |
//creat force Subtraction policy |
31 |
+ |
data = info->getProperty(ZCONSFORCEPOLICY_ID); |
32 |
+ |
if(!data){ |
33 |
|
sprintf( painCave.errMsg, |
34 |
< |
"ZConstraint error: If you use an ZConstraint\n" |
35 |
< |
" , you must set index of z-constraint molecules.\n"); |
36 |
< |
painCave.isFatal = 1; |
37 |
< |
simError(); |
34 |
> |
"ZConstraint Warning: User does not set force Subtraction policy, " |
35 |
> |
"PolicyByMass is used\n"); |
36 |
> |
painCave.isFatal = 0; |
37 |
> |
simError(); |
38 |
> |
|
39 |
> |
forcePolicy = (ForceSubtractionPolicy*) new PolicyByMass(this); |
40 |
|
} |
41 |
|
else{ |
42 |
< |
index = dynamic_cast<IndexData*>(data); |
42 |
> |
policy = dynamic_cast<StringData*>(data); |
43 |
|
|
44 |
< |
if(!index){ |
28 |
< |
|
44 |
> |
if(!policy){ |
45 |
|
sprintf( painCave.errMsg, |
46 |
< |
"ZConstraint error: Can not get property from SimInfo\n"); |
47 |
< |
painCave.isFatal = 1; |
48 |
< |
simError(); |
49 |
< |
|
46 |
> |
"ZConstraint Error: Convertion from GenericData to StringData failure, " |
47 |
> |
"PolicyByMass is used\n"); |
48 |
> |
painCave.isFatal = 0; |
49 |
> |
simError(); |
50 |
> |
|
51 |
> |
forcePolicy = (ForceSubtractionPolicy*) new PolicyByMass(this); |
52 |
|
} |
53 |
|
else{ |
54 |
< |
|
55 |
< |
indexOfAllZConsMols = index->getIndexData(); |
56 |
< |
|
57 |
< |
//the maximum value of index is the last one(we sorted the index data in SimSetup.cpp) |
58 |
< |
int maxIndex; |
41 |
< |
int totalNumMol; |
42 |
< |
|
43 |
< |
maxIndex = indexOfAllZConsMols[indexOfAllZConsMols.size() - 1]; |
44 |
< |
|
45 |
< |
#ifndef IS_MPI |
46 |
< |
totalNumMol = nMols; |
47 |
< |
#else |
48 |
< |
totalNumMol = mpiSim->getTotNmol(); |
49 |
< |
#endif |
50 |
< |
|
51 |
< |
if(maxIndex > totalNumMol - 1){ |
54 |
> |
if(policy->getData() == "BYNUMBER") |
55 |
> |
forcePolicy = (ForceSubtractionPolicy*) new PolicyByNumber(this); |
56 |
> |
else if(policy->getData() == "BYMASS") |
57 |
> |
forcePolicy = (ForceSubtractionPolicy*) new PolicyByMass(this); |
58 |
> |
else{ |
59 |
|
sprintf( painCave.errMsg, |
60 |
< |
"ZConstraint error: index is out of range\n"); |
61 |
< |
painCave.isFatal = 1; |
62 |
< |
simError(); |
63 |
< |
|
64 |
< |
} |
65 |
< |
|
60 |
> |
"ZConstraint Warning: unknown force Subtraction policy, " |
61 |
> |
"PolicyByMass is used\n"); |
62 |
> |
painCave.isFatal = 0; |
63 |
> |
simError(); |
64 |
> |
forcePolicy = (ForceSubtractionPolicy*) new PolicyByMass(this); |
65 |
> |
} |
66 |
|
} |
60 |
– |
|
67 |
|
} |
68 |
|
|
63 |
– |
//retrive sample time of z-contraint |
64 |
– |
data = info->getProperty("zconstime"); |
69 |
|
|
70 |
+ |
//retrieve sample time of z-contraint |
71 |
+ |
data = info->getProperty(ZCONSTIME_ID); |
72 |
+ |
|
73 |
|
if(!data) { |
74 |
|
|
75 |
|
sprintf( painCave.errMsg, |
96 |
|
|
97 |
|
} |
98 |
|
|
99 |
< |
|
100 |
< |
//retrive output filename of z force |
94 |
< |
data = info->getProperty("zconsfilename"); |
99 |
> |
//retrieve output filename of z force |
100 |
> |
data = info->getProperty(ZCONSFILENAME_ID); |
101 |
|
if(!data) { |
102 |
|
|
103 |
|
|
110 |
|
} |
111 |
|
else{ |
112 |
|
|
113 |
< |
filename = dynamic_cast<StringData*>(data); |
113 |
> |
filename = dynamic_cast<StringData*>(data); |
114 |
|
|
115 |
|
if(!filename){ |
116 |
|
|
124 |
|
this->zconsOutput = filename->getData(); |
125 |
|
} |
126 |
|
|
121 |
– |
|
127 |
|
} |
128 |
|
|
129 |
< |
|
130 |
< |
//calculate reference z coordinate for z-constraint molecules |
126 |
< |
double totalMass_local; |
127 |
< |
double totalMass; |
128 |
< |
double totalMZ_local; |
129 |
< |
double totalMZ; |
130 |
< |
double massOfUncons_local; |
131 |
< |
double massOfCurMol; |
132 |
< |
double COM[3]; |
129 |
> |
//retrieve tolerance for z-constraint molecuels |
130 |
> |
data = info->getProperty(ZCONSTOL_ID); |
131 |
|
|
132 |
< |
totalMass_local = 0; |
133 |
< |
totalMass = 0; |
134 |
< |
totalMZ_local = 0; |
135 |
< |
totalMZ = 0; |
136 |
< |
massOfUncons_local = 0; |
137 |
< |
|
140 |
< |
|
141 |
< |
for(int i = 0; i < nMols; i++){ |
142 |
< |
massOfCurMol = molecules[i].getTotalMass(); |
143 |
< |
molecules[i].getCOM(COM); |
144 |
< |
|
145 |
< |
totalMass_local += massOfCurMol; |
146 |
< |
totalMZ_local += massOfCurMol * COM[2]; |
147 |
< |
|
148 |
< |
if(isZConstraintMol(&molecules[i]) == -1){ |
149 |
< |
|
150 |
< |
massOfUncons_local += massOfCurMol; |
151 |
< |
} |
152 |
< |
|
132 |
> |
if(!data) { |
133 |
> |
|
134 |
> |
sprintf( painCave.errMsg, |
135 |
> |
"ZConstraint error: can not get tolerance \n"); |
136 |
> |
painCave.isFatal = 1; |
137 |
> |
simError(); |
138 |
|
} |
139 |
+ |
else{ |
140 |
|
|
141 |
< |
|
142 |
< |
#ifdef IS_MPI |
143 |
< |
MPI_Allreduce(&totalMass_local, &totalMass, 1, MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
158 |
< |
MPI_Allreduce(&totalMZ_local, &totalMZ, 1, MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
159 |
< |
MPI_Allreduce(&massOfUncons_local, &totalMassOfUncons, 1, MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
160 |
< |
#else |
161 |
< |
totalMass = totalMass_local; |
162 |
< |
totalMZ = totalMZ_local; |
163 |
< |
totalMassOfUncons = massOfUncons_local; |
164 |
< |
#endif |
141 |
> |
tolerance = dynamic_cast<DoubleData*>(data); |
142 |
> |
|
143 |
> |
if(!tolerance){ |
144 |
|
|
145 |
< |
double zsys; |
146 |
< |
zsys = totalMZ / totalMass; |
147 |
< |
|
148 |
< |
#ifndef IS_MPI |
149 |
< |
for(int i = 0; i < nMols; i++){ |
171 |
< |
|
172 |
< |
if(isZConstraintMol(&molecules[i]) > -1 ){ |
173 |
< |
molecules[i].getCOM(COM); |
174 |
< |
allRefZ.push_back(COM[2] - zsys); |
145 |
> |
sprintf( painCave.errMsg, |
146 |
> |
"ZConstraint error: Can not get property from SimInfo\n"); |
147 |
> |
painCave.isFatal = 1; |
148 |
> |
simError(); |
149 |
> |
|
150 |
|
} |
151 |
< |
|
151 |
> |
else{ |
152 |
> |
this->zconsTol = tolerance->getData(); |
153 |
> |
} |
154 |
> |
|
155 |
|
} |
178 |
– |
#else |
179 |
– |
|
180 |
– |
int whichNode; |
181 |
– |
enum CommType { RequestMolZPos, EndOfRequest} status; |
182 |
– |
//int status; |
183 |
– |
double zpos; |
184 |
– |
int localIndex; |
185 |
– |
MPI_Status ierr; |
186 |
– |
int tag = 0; |
156 |
|
|
157 |
< |
if(worldRank == 0){ |
158 |
< |
|
159 |
< |
int globalIndexOfCurMol; |
160 |
< |
int *MolToProcMap; |
161 |
< |
MolToProcMap = mpiSim->getMolToProcMap(); |
162 |
< |
|
163 |
< |
for(int i = 0; i < indexOfAllZConsMols.size(); i++){ |
164 |
< |
|
165 |
< |
whichNode = MolToProcMap[indexOfAllZConsMols[i]]; |
197 |
< |
globalIndexOfCurMol = indexOfAllZConsMols[i]; |
198 |
< |
|
199 |
< |
if(whichNode == 0){ |
200 |
< |
|
201 |
< |
for(int j = 0; j < nMols; j++) |
202 |
< |
if(molecules[j].getGlobalIndex() == globalIndexOfCurMol){ |
203 |
< |
localIndex = j; |
204 |
< |
break; |
205 |
< |
} |
206 |
< |
|
207 |
< |
molecules[localIndex].getCOM(COM); |
208 |
< |
allRefZ.push_back(COM[2] - zsys); |
209 |
< |
|
210 |
< |
} |
211 |
< |
else{ |
212 |
< |
status = RequestMolZPos; |
213 |
< |
MPI_Send(&status, 1, MPI_INT, whichNode, tag, MPI_COMM_WORLD); |
214 |
< |
MPI_Send(&globalIndexOfCurMol, 1, MPI_INT, whichNode, tag, MPI_COMM_WORLD); |
215 |
< |
MPI_Recv(&zpos, 1, MPI_DOUBLE_PRECISION, whichNode, tag, MPI_COMM_WORLD, &ierr); |
216 |
< |
|
217 |
< |
allRefZ.push_back(zpos - zsys); |
218 |
< |
|
219 |
< |
} |
220 |
< |
|
221 |
< |
} //End of Request Loop |
222 |
< |
|
223 |
< |
//Send ending request message to slave nodes |
224 |
< |
status = EndOfRequest; |
225 |
< |
for(int i =1; i < mpiSim->getNumberProcessors(); i++) |
226 |
< |
MPI_Send(&status, 1, MPI_INT, i, tag, MPI_COMM_WORLD); |
227 |
< |
|
157 |
> |
//retrieve index of z-constraint molecules |
158 |
> |
data = info->getProperty(ZCONSPARADATA_ID); |
159 |
> |
if(!data) { |
160 |
> |
|
161 |
> |
sprintf( painCave.errMsg, |
162 |
> |
"ZConstraint error: If you use an ZConstraint\n" |
163 |
> |
" , you must set index of z-constraint molecules.\n"); |
164 |
> |
painCave.isFatal = 1; |
165 |
> |
simError(); |
166 |
|
} |
167 |
|
else{ |
168 |
|
|
169 |
< |
int whichMol; |
170 |
< |
bool done = false; |
169 |
> |
zConsParaData = dynamic_cast<ZConsParaData*>(data); |
170 |
> |
|
171 |
> |
if(!zConsParaData){ |
172 |
|
|
173 |
< |
while (!done){ |
174 |
< |
|
175 |
< |
MPI_Recv(&status, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &ierr); |
173 |
> |
sprintf( painCave.errMsg, |
174 |
> |
"ZConstraint error: Can not get parameters of zconstraint method from SimInfo\n"); |
175 |
> |
painCave.isFatal = 1; |
176 |
> |
simError(); |
177 |
|
|
238 |
– |
switch (status){ |
239 |
– |
|
240 |
– |
case RequestMolZPos : |
241 |
– |
|
242 |
– |
MPI_Recv(&whichMol, 1, MPI_INT, 0, tag, MPI_COMM_WORLD,&ierr); |
243 |
– |
|
244 |
– |
for(int i = 0; i < nMols; i++) |
245 |
– |
if(molecules[i].getGlobalIndex() == whichMol){ |
246 |
– |
localIndex = i; |
247 |
– |
break; |
248 |
– |
} |
249 |
– |
|
250 |
– |
molecules[localIndex].getCOM(COM); |
251 |
– |
zpos = COM[2]; |
252 |
– |
MPI_Send(&zpos, 1, MPI_DOUBLE_PRECISION, 0, tag, MPI_COMM_WORLD); |
253 |
– |
|
254 |
– |
break; |
255 |
– |
|
256 |
– |
case EndOfRequest : |
257 |
– |
|
258 |
– |
done = true; |
259 |
– |
break; |
260 |
– |
} |
261 |
– |
|
178 |
|
} |
179 |
< |
|
180 |
< |
} |
179 |
> |
else{ |
180 |
> |
|
181 |
> |
parameters = zConsParaData->getData(); |
182 |
|
|
183 |
< |
//Brocast the allRefZ to slave nodes; |
184 |
< |
double* allRefZBuf; |
185 |
< |
int nZConsMols; |
269 |
< |
nZConsMols = indexOfAllZConsMols.size(); |
270 |
< |
|
271 |
< |
allRefZBuf = new double[nZConsMols]; |
272 |
< |
|
273 |
< |
if(worldRank == 0){ |
183 |
> |
//check the range of zconsIndex |
184 |
> |
//and the minimum value of index is the first one (we already sorted the data) |
185 |
> |
//the maximum value of index is the last one |
186 |
|
|
187 |
< |
for(int i = 0; i < nZConsMols; i++) |
188 |
< |
allRefZBuf[i] = allRefZ[i]; |
189 |
< |
} |
190 |
< |
|
191 |
< |
MPI_Bcast(allRefZBuf, nZConsMols, MPI_DOUBLE_PRECISION, 0, MPI_COMM_WORLD); |
192 |
< |
|
193 |
< |
if(worldRank != 0){ |
194 |
< |
|
195 |
< |
for(int i = 0; i < nZConsMols; i++) |
196 |
< |
allRefZ.push_back(allRefZBuf[i]); |
197 |
< |
} |
198 |
< |
|
199 |
< |
delete[] allRefZBuf; |
187 |
> |
int maxIndex; |
188 |
> |
int minIndex; |
189 |
> |
int totalNumMol; |
190 |
> |
|
191 |
> |
minIndex = (*parameters)[0].zconsIndex; |
192 |
> |
if(minIndex < 0){ |
193 |
> |
sprintf( painCave.errMsg, |
194 |
> |
"ZConstraint error: index is out of range\n"); |
195 |
> |
painCave.isFatal = 1; |
196 |
> |
simError(); |
197 |
> |
} |
198 |
> |
|
199 |
> |
maxIndex = (*parameters)[parameters->size() - 1].zconsIndex; |
200 |
> |
|
201 |
> |
#ifndef IS_MPI |
202 |
> |
totalNumMol = nMols; |
203 |
> |
#else |
204 |
> |
totalNumMol = mpiSim->getTotNmol(); |
205 |
> |
#endif |
206 |
> |
|
207 |
> |
if(maxIndex > totalNumMol - 1){ |
208 |
> |
sprintf( painCave.errMsg, |
209 |
> |
"ZConstraint error: index is out of range\n"); |
210 |
> |
painCave.isFatal = 1; |
211 |
> |
simError(); |
212 |
> |
} |
213 |
> |
|
214 |
> |
//if user does not specify the zpos for the zconstraint molecule |
215 |
> |
//its initial z coordinate will be used as default |
216 |
> |
for(int i = 0; i < parameters->size(); i++){ |
217 |
> |
|
218 |
> |
if(!(*parameters)[i].havingZPos){ |
219 |
> |
#ifndef IS_MPI |
220 |
> |
for(int j = 0; j < nMols; j++){ |
221 |
> |
if (molecules[j].getGlobalIndex() == (*parameters)[i].zconsIndex){ |
222 |
> |
molecules[j].getCOM(COM); |
223 |
> |
break; |
224 |
> |
} |
225 |
> |
} |
226 |
> |
#else |
227 |
> |
//query which processor current zconstraint molecule belongs to |
228 |
> |
int *MolToProcMap; |
229 |
> |
int whichNode; |
230 |
> |
double initZPos; |
231 |
> |
MolToProcMap = mpiSim->getMolToProcMap(); |
232 |
> |
whichNode = MolToProcMap[(*parameters)[i].zconsIndex]; |
233 |
> |
|
234 |
> |
//broadcast the zpos of current z-contraint molecule |
235 |
> |
//the node which contain this |
236 |
> |
|
237 |
> |
if (worldRank == whichNode ){ |
238 |
> |
|
239 |
> |
for(int j = 0; j < nMols; j++) |
240 |
> |
if (molecules[j].getGlobalIndex() == (*parameters)[i].zconsIndex){ |
241 |
> |
molecules[j].getCOM(COM); |
242 |
> |
break; |
243 |
> |
} |
244 |
> |
|
245 |
> |
} |
246 |
> |
|
247 |
> |
MPI_Bcast(&COM[whichDirection], 1, MPI_DOUBLE_PRECISION, whichNode, MPI_COMM_WORLD); |
248 |
|
#endif |
249 |
+ |
|
250 |
+ |
(*parameters)[i].zPos = COM[whichDirection]; |
251 |
|
|
252 |
< |
|
252 |
> |
sprintf( painCave.errMsg, |
253 |
> |
"ZConstraint warning: Does not specify zpos for z-constraint molecule " |
254 |
> |
"initial z coornidate will be used \n"); |
255 |
> |
painCave.isFatal = 0; |
256 |
> |
simError(); |
257 |
> |
|
258 |
> |
} |
259 |
> |
} |
260 |
> |
|
261 |
> |
}//end if (!zConsParaData) |
262 |
> |
}//end if (!data) |
263 |
> |
|
264 |
> |
// |
265 |
|
#ifdef IS_MPI |
266 |
|
update(); |
267 |
|
#else |
268 |
|
int searchResult; |
269 |
< |
|
296 |
< |
refZ = allRefZ; |
297 |
< |
|
269 |
> |
|
270 |
|
for(int i = 0; i < nMols; i++){ |
271 |
|
|
272 |
|
searchResult = isZConstraintMol(&molecules[i]); |
275 |
|
|
276 |
|
zconsMols.push_back(&molecules[i]); |
277 |
|
massOfZConsMols.push_back(molecules[i].getTotalMass()); |
278 |
< |
|
278 |
> |
|
279 |
> |
zPos.push_back((*parameters)[searchResult].zPos); |
280 |
> |
// cout << "index: "<< (*parameters)[searchResult].zconsIndex |
281 |
> |
// <<"\tzPos = " << (*parameters)[searchResult].zPos << endl; |
282 |
> |
|
283 |
> |
kz.push_back((*parameters)[searchResult]. kRatio * zForceConst); |
284 |
|
molecules[i].getCOM(COM); |
285 |
|
} |
286 |
|
else |
293 |
|
} |
294 |
|
|
295 |
|
fz = new double[zconsMols.size()]; |
296 |
+ |
curZPos = new double[zconsMols.size()]; |
297 |
|
indexOfZConsMols = new int [zconsMols.size()]; |
298 |
|
|
299 |
< |
if(!fz || !indexOfZConsMols){ |
299 |
> |
if(!fz || !curZPos || !indexOfZConsMols){ |
300 |
|
sprintf( painCave.errMsg, |
301 |
|
"Memory allocation failure in class Zconstraint\n"); |
302 |
|
painCave.isFatal = 1; |
303 |
|
simError(); |
304 |
|
} |
305 |
|
|
306 |
< |
for(int i = 0; i < zconsMols.size(); i++) |
306 |
> |
//determine the states of z-constraint molecules |
307 |
> |
for(int i = 0; i < zconsMols.size(); i++){ |
308 |
|
indexOfZConsMols[i] = zconsMols[i]->getGlobalIndex(); |
309 |
+ |
|
310 |
+ |
zconsMols[i]->getCOM(COM); |
311 |
+ |
if (fabs(zPos[i] - COM[whichDirection]) < zconsTol) |
312 |
+ |
states.push_back(zcsFixed); |
313 |
+ |
else |
314 |
+ |
states.push_back(zcsMoving); |
315 |
+ |
} |
316 |
|
|
317 |
|
#endif |
318 |
+ |
|
319 |
+ |
//get total masss of unconstraint molecules |
320 |
+ |
double totalMassOfUncons_local; |
321 |
+ |
totalMassOfUncons_local = 0; |
322 |
|
|
323 |
< |
fzOut = new ZConsWriter(zconsOutput.c_str()); |
324 |
< |
|
325 |
< |
if(!fzOut){ |
326 |
< |
sprintf( painCave.errMsg, |
327 |
< |
"Memory allocation failure in class Zconstraint\n"); |
328 |
< |
painCave.isFatal = 1; |
329 |
< |
simError(); |
330 |
< |
} |
331 |
< |
|
332 |
< |
fzOut->writeRefZ(indexOfAllZConsMols, allRefZ); |
323 |
> |
for(int i = 0; i < unconsMols.size(); i++) |
324 |
> |
totalMassOfUncons_local += unconsMols[i]->getTotalMass(); |
325 |
> |
|
326 |
> |
#ifndef IS_MPI |
327 |
> |
totalMassOfUncons = totalMassOfUncons_local; |
328 |
> |
#else |
329 |
> |
MPI_Allreduce(&totalMassOfUncons_local, &totalMassOfUncons, 1, |
330 |
> |
MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD); |
331 |
> |
#endif |
332 |
> |
|
333 |
> |
//get total number of unconstrained atoms |
334 |
> |
int nUnconsAtoms_local; |
335 |
> |
nUnconsAtoms_local = 0; |
336 |
> |
for(int i = 0; i < unconsMols.size(); i++) |
337 |
> |
nUnconsAtoms_local += unconsMols[i]->getNAtoms(); |
338 |
> |
|
339 |
> |
#ifndef IS_MPI |
340 |
> |
totNumOfUnconsAtoms = nUnconsAtoms_local; |
341 |
> |
#else |
342 |
> |
MPI_Allreduce(&nUnconsAtoms_local, &totNumOfUnconsAtoms, 1, |
343 |
> |
MPI_INT,MPI_SUM, MPI_COMM_WORLD); |
344 |
> |
#endif |
345 |
> |
|
346 |
> |
forcePolicy->update(); |
347 |
|
} |
348 |
|
|
349 |
|
template<typename T> ZConstraint<T>::~ZConstraint() |
350 |
|
{ |
351 |
|
if(fz) |
352 |
|
delete[] fz; |
353 |
+ |
|
354 |
+ |
if(curZPos) |
355 |
+ |
delete[] curZPos; |
356 |
|
|
357 |
|
if(indexOfZConsMols) |
358 |
|
delete[] indexOfZConsMols; |
359 |
|
|
360 |
|
if(fzOut) |
361 |
|
delete fzOut; |
362 |
+ |
|
363 |
+ |
if(forcePolicy) |
364 |
+ |
delete forcePolicy; |
365 |
|
} |
366 |
|
|
367 |
+ |
|
368 |
+ |
/** |
369 |
+ |
* |
370 |
+ |
*/ |
371 |
+ |
|
372 |
|
#ifdef IS_MPI |
373 |
|
template<typename T> void ZConstraint<T>::update() |
374 |
|
{ |
377 |
|
|
378 |
|
zconsMols.clear(); |
379 |
|
massOfZConsMols.clear(); |
380 |
< |
refZ.clear(); |
380 |
> |
zPos.clear(); |
381 |
> |
kz.clear(); |
382 |
|
|
383 |
|
unconsMols.clear(); |
384 |
|
massOfUnconsMols.clear(); |
392 |
|
if(index > -1){ |
393 |
|
|
394 |
|
zconsMols.push_back(&molecules[i]); |
395 |
+ |
zPos.push_back((*parameters)[index].zPos); |
396 |
+ |
kz.push_back((*parameters)[index].kRatio * zForceConst); |
397 |
|
massOfZConsMols.push_back(molecules[i].getTotalMass()); |
398 |
|
|
399 |
|
molecules[i].getCOM(COM); |
382 |
– |
refZ.push_back(allRefZ[index]); |
400 |
|
} |
401 |
|
else |
402 |
|
{ |
406 |
|
|
407 |
|
} |
408 |
|
} |
409 |
+ |
|
410 |
+ |
//determine the states of z-constraint molecules |
411 |
+ |
for(int i = 0; i < zconsMols.size(); i++){ |
412 |
+ |
zconsMols[i]->getCOM(COM); |
413 |
+ |
if (fabs(zPos[i] - COM[whichDirection]) < zconsTol) |
414 |
+ |
states.push_back(zcsFixed); |
415 |
+ |
else |
416 |
+ |
states.push_back(zcsMoving); |
417 |
+ |
} |
418 |
+ |
|
419 |
|
|
420 |
|
//The reason to declare fz and indexOfZconsMols as pointer to array is |
421 |
|
// that we want to make the MPI communication simple |
422 |
|
if(fz) |
423 |
|
delete[] fz; |
424 |
+ |
|
425 |
+ |
if(curZPos) |
426 |
+ |
delete[] curZPos; |
427 |
|
|
428 |
|
if(indexOfZConsMols) |
429 |
|
delete[] indexOfZConsMols; |
430 |
|
|
431 |
|
if (zconsMols.size() > 0){ |
432 |
|
fz = new double[zconsMols.size()]; |
433 |
+ |
curZPos = new double[zconsMols.size()]; |
434 |
|
indexOfZConsMols = new int[zconsMols.size()]; |
435 |
|
|
436 |
< |
if(!fz || !indexOfZConsMols){ |
436 |
> |
if(!fz || !curZPos || !indexOfZConsMols){ |
437 |
|
sprintf( painCave.errMsg, |
438 |
|
"Memory allocation failure in class Zconstraint\n"); |
439 |
|
painCave.isFatal = 1; |
447 |
|
} |
448 |
|
else{ |
449 |
|
fz = NULL; |
450 |
+ |
curZPos = NULL; |
451 |
|
indexOfZConsMols = NULL; |
452 |
|
} |
453 |
|
|
454 |
+ |
// |
455 |
+ |
forcePolicy->update(); |
456 |
+ |
|
457 |
|
} |
458 |
|
|
459 |
|
#endif |
460 |
|
|
461 |
< |
/** Function Name: isZConstraintMol |
462 |
< |
** Parameter |
463 |
< |
** Molecule* mol |
464 |
< |
** Return value: |
465 |
< |
** -1, if the molecule is not z-constraint molecule, |
466 |
< |
** other non-negative values, its index in indexOfAllZConsMols vector |
461 |
> |
/** |
462 |
> |
* Function Name: isZConstraintMol |
463 |
> |
* Parameter |
464 |
> |
* Molecule* mol |
465 |
> |
* Return value: |
466 |
> |
* -1, if the molecule is not z-constraint molecule, |
467 |
> |
* other non-negative values, its index in indexOfAllZConsMols vector |
468 |
|
*/ |
469 |
|
|
470 |
|
template<typename T> int ZConstraint<T>::isZConstraintMol(Molecule* mol) |
477 |
|
index = mol->getGlobalIndex(); |
478 |
|
|
479 |
|
low = 0; |
480 |
< |
high = indexOfAllZConsMols.size() - 1; |
480 |
> |
high = parameters->size() - 1; |
481 |
|
|
482 |
|
//Binary Search (we have sorted the array) |
483 |
|
while(low <= high){ |
484 |
|
mid = (low + high) /2; |
485 |
< |
if (indexOfAllZConsMols[mid] == index) |
485 |
> |
if ((*parameters)[mid].zconsIndex == index) |
486 |
|
return mid; |
487 |
< |
else if (indexOfAllZConsMols[mid] > index ) |
487 |
> |
else if ((*parameters)[mid].zconsIndex > index ) |
488 |
|
high = mid -1; |
489 |
|
else |
490 |
|
low = mid + 1; |
493 |
|
return -1; |
494 |
|
} |
495 |
|
|
496 |
< |
/** Function Name: integrateStep |
497 |
< |
** Parameter: |
498 |
< |
** int calcPot; |
499 |
< |
** int calcStress; |
464 |
< |
** Description: |
465 |
< |
** Advance One Step. |
466 |
< |
** Memo: |
467 |
< |
** The best way to implement z-constraint is to override integrateStep |
468 |
< |
** Overriding constrainB is not a good choice, since in integrateStep, |
469 |
< |
** constrainB is invoked by below line, |
470 |
< |
** if(nConstrained) constrainB(); |
471 |
< |
** For instance, we would like to apply z-constraint without bond contrain, |
472 |
< |
** In that case, if we override constrainB, Z-constrain method will never be executed; |
473 |
< |
*/ |
474 |
< |
template<typename T> void ZConstraint<T>::integrateStep( int calcPot, int calcStress ) |
475 |
< |
{ |
476 |
< |
T::integrateStep( calcPot, calcStress ); |
477 |
< |
resetZ(); |
496 |
> |
template<typename T> void ZConstraint<T>::integrate(){ |
497 |
> |
|
498 |
> |
// creat zconsWriter |
499 |
> |
fzOut = new ZConsWriter(zconsOutput.c_str(), parameters); |
500 |
|
|
501 |
< |
double currZConsTime = 0; |
502 |
< |
|
503 |
< |
//write out forces of z constraint |
504 |
< |
if( info->getTime() >= currZConsTime){ |
505 |
< |
fzOut->writeFZ(info->getTime(), zconsMols.size(),indexOfZConsMols, fz); |
506 |
< |
} |
507 |
< |
} |
501 |
> |
if(!fzOut){ |
502 |
> |
sprintf( painCave.errMsg, |
503 |
> |
"Memory allocation failure in class Zconstraint\n"); |
504 |
> |
painCave.isFatal = 1; |
505 |
> |
simError(); |
506 |
> |
} |
507 |
> |
|
508 |
> |
//zero out the velocities of center of mass of unconstrained molecules |
509 |
> |
//and the velocities of center of mass of every single z-constrained molecueles |
510 |
> |
zeroOutVel(); |
511 |
|
|
512 |
< |
/** Function Name: resetZ |
513 |
< |
** Description: |
514 |
< |
** Reset the z coordinates |
490 |
< |
*/ |
512 |
> |
curZconsTime = zconsTime + info->getTime(); |
513 |
> |
|
514 |
> |
T::integrate(); |
515 |
|
|
516 |
< |
template<typename T> void ZConstraint<T>::resetZ() |
517 |
< |
{ |
518 |
< |
double deltaZ; |
519 |
< |
double mzOfZCons; //total sum of m*z of z-constrain molecules |
520 |
< |
double mzOfUncons; //total sum of m*z of unconstrain molecuels; |
521 |
< |
double totalMZOfZCons; |
522 |
< |
double totalMZOfUncons; |
523 |
< |
double COM[3]; |
516 |
> |
} |
517 |
> |
|
518 |
> |
|
519 |
> |
/** |
520 |
> |
* |
521 |
> |
* |
522 |
> |
* |
523 |
> |
* |
524 |
> |
*/ |
525 |
> |
template<typename T> void ZConstraint<T>::calcForce(int calcPot, int calcStress){ |
526 |
|
double zsys; |
527 |
< |
Atom** zconsAtoms; |
527 |
> |
double COM[3]; |
528 |
> |
double force[3]; |
529 |
> |
double zSysCOMVel; |
530 |
|
|
531 |
< |
mzOfZCons = 0; |
532 |
< |
mzOfUncons = 0; |
531 |
> |
T::calcForce(calcPot, calcStress); |
532 |
> |
|
533 |
> |
if (checkZConsState()){ |
534 |
> |
zeroOutVel(); |
535 |
> |
forcePolicy->update(); |
536 |
> |
} |
537 |
|
|
538 |
< |
for(int i = 0; i < zconsMols.size(); i++){ |
539 |
< |
mzOfZCons += massOfZConsMols[i] * refZ[i]; |
540 |
< |
} |
538 |
> |
zsys = calcZSys(); |
539 |
> |
zSysCOMVel = calcSysCOMVel(); |
540 |
> |
#ifdef IS_MPI |
541 |
> |
if(worldRank == 0){ |
542 |
> |
#endif |
543 |
> |
//cout << "---------------------------------------------------------------------" <<endl; |
544 |
> |
//cout << "current time: " << info->getTime() << endl; |
545 |
> |
//cout << "center of mass at z: " << zsys << endl; |
546 |
> |
//cout << "before calcForce, the COMVel of system is " << zSysCOMVel <<endl; |
547 |
|
|
548 |
|
#ifdef IS_MPI |
549 |
< |
MPI_Allreduce(&mzOfZCons, &totalMZOfZCons, 1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
512 |
< |
#else |
513 |
< |
totalMZOfZCons = mzOfZCons; |
549 |
> |
} |
550 |
|
#endif |
551 |
|
|
552 |
< |
for(int i = 0; i < unconsMols.size(); i++){ |
553 |
< |
unconsMols[i]->getCOM(COM); |
554 |
< |
mzOfUncons += massOfUnconsMols[i] * COM[2]; |
552 |
> |
//do zconstraint force; |
553 |
> |
if (haveFixedZMols()) |
554 |
> |
this->doZconstraintForce(); |
555 |
> |
|
556 |
> |
//use harmonical poteintial to move the molecules to the specified positions |
557 |
> |
if (haveMovingZMols()) |
558 |
> |
this->doHarmonic(); |
559 |
> |
|
560 |
> |
//write out forces and current positions of z-constraint molecules |
561 |
> |
if(info->getTime() >= curZconsTime){ |
562 |
> |
for(int i = 0; i < zconsMols.size(); i++){ |
563 |
> |
zconsMols[i]->getCOM(COM); |
564 |
> |
curZPos[i] = COM[whichDirection]; |
565 |
> |
|
566 |
> |
//if the z-constraint molecule is still moving, just record its force |
567 |
> |
if(states[i] == zcsMoving){ |
568 |
> |
fz[i] = 0; |
569 |
> |
Atom** movingZAtoms; |
570 |
> |
movingZAtoms = zconsMols[i]->getMyAtoms(); |
571 |
> |
for(int j = 0; j < zconsMols[i]->getNAtoms(); j++){ |
572 |
> |
movingZAtoms[j]->getFrc(force); |
573 |
> |
fz[i] += force[whichDirection]; |
574 |
> |
} |
575 |
> |
} |
576 |
> |
} |
577 |
> |
fzOut->writeFZ(info->getTime(), zconsMols.size(), indexOfZConsMols, fz, curZPos); |
578 |
> |
curZconsTime += zconsTime; |
579 |
|
} |
580 |
< |
|
580 |
> |
|
581 |
> |
zSysCOMVel = calcSysCOMVel(); |
582 |
|
#ifdef IS_MPI |
583 |
< |
MPI_Allreduce(&mzOfUncons, &totalMZOfUncons, 1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
584 |
< |
#else |
585 |
< |
totalMZOfUncons = mzOfUncons; |
586 |
< |
#endif |
587 |
< |
|
588 |
< |
zsys = (totalMZOfZCons + totalMZOfUncons) /totalMassOfUncons; |
583 |
> |
if(worldRank == 0){ |
584 |
> |
#endif |
585 |
> |
//cout << "after calcForce, the COMVel of system is " << zSysCOMVel <<endl; |
586 |
> |
#ifdef IS_MPI |
587 |
> |
} |
588 |
> |
#endif |
589 |
> |
} |
590 |
> |
|
591 |
> |
|
592 |
> |
/** |
593 |
> |
* |
594 |
> |
*/ |
595 |
|
|
596 |
< |
for(int i = 0; i < zconsMols.size(); i++){ |
596 |
> |
template<typename T> double ZConstraint<T>::calcZSys() |
597 |
> |
{ |
598 |
> |
//calculate reference z coordinate for z-constraint molecules |
599 |
> |
double totalMass_local; |
600 |
> |
double totalMass; |
601 |
> |
double totalMZ_local; |
602 |
> |
double totalMZ; |
603 |
> |
double massOfCurMol; |
604 |
> |
double COM[3]; |
605 |
> |
|
606 |
> |
totalMass_local = 0; |
607 |
> |
totalMZ_local = 0; |
608 |
> |
|
609 |
> |
for(int i = 0; i < nMols; i++){ |
610 |
> |
massOfCurMol = molecules[i].getTotalMass(); |
611 |
> |
molecules[i].getCOM(COM); |
612 |
> |
|
613 |
> |
totalMass_local += massOfCurMol; |
614 |
> |
totalMZ_local += massOfCurMol * COM[whichDirection]; |
615 |
> |
|
616 |
> |
} |
617 |
> |
|
618 |
> |
|
619 |
> |
#ifdef IS_MPI |
620 |
> |
MPI_Allreduce(&totalMass_local, &totalMass, 1, |
621 |
> |
MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD); |
622 |
> |
MPI_Allreduce(&totalMZ_local, &totalMZ, 1, |
623 |
> |
MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD); |
624 |
> |
#else |
625 |
> |
totalMass = totalMass_local; |
626 |
> |
totalMZ = totalMZ_local; |
627 |
> |
#endif |
628 |
> |
|
629 |
> |
double zsys; |
630 |
> |
zsys = totalMZ / totalMass; |
631 |
> |
|
632 |
> |
return zsys; |
633 |
> |
} |
634 |
> |
|
635 |
> |
/** |
636 |
> |
* |
637 |
> |
*/ |
638 |
> |
template<typename T> void ZConstraint<T>::thermalize( void ){ |
639 |
> |
|
640 |
> |
T::thermalize(); |
641 |
> |
zeroOutVel(); |
642 |
> |
} |
643 |
> |
|
644 |
> |
/** |
645 |
> |
* |
646 |
> |
*/ |
647 |
> |
|
648 |
> |
template<typename T> void ZConstraint<T>::zeroOutVel(){ |
649 |
> |
|
650 |
> |
Atom** fixedZAtoms; |
651 |
> |
double COMvel[3]; |
652 |
> |
double vel[3]; |
653 |
> |
double zSysCOMVel; |
654 |
> |
|
655 |
> |
//zero out the velocities of center of mass of fixed z-constrained molecules |
656 |
> |
|
657 |
> |
for(int i = 0; i < zconsMols.size(); i++){ |
658 |
> |
|
659 |
> |
if (states[i] == zcsFixed){ |
660 |
> |
|
661 |
> |
zconsMols[i]->getCOMvel(COMvel); |
662 |
> |
//cout << "before resetting " << indexOfZConsMols[i] <<"'s vz is " << COMvel[whichDirection] << endl; |
663 |
> |
|
664 |
> |
fixedZAtoms = zconsMols[i]->getMyAtoms(); |
665 |
> |
|
666 |
> |
for(int j =0; j < zconsMols[i]->getNAtoms(); j++){ |
667 |
> |
fixedZAtoms[j]->getVel(vel); |
668 |
> |
vel[whichDirection] -= COMvel[whichDirection]; |
669 |
> |
fixedZAtoms[j]->setVel(vel); |
670 |
> |
} |
671 |
> |
|
672 |
> |
zconsMols[i]->getCOMvel(COMvel); |
673 |
> |
//cout << "after resetting " << indexOfZConsMols[i] <<"'s vz is " << COMvel[whichDirection] << endl; |
674 |
> |
} |
675 |
> |
|
676 |
> |
} |
677 |
> |
|
678 |
> |
//cout << "before resetting the COMVel of moving molecules is " << calcMovingMolsCOMVel() <<endl; |
679 |
> |
|
680 |
> |
zSysCOMVel = calcSysCOMVel(); |
681 |
> |
#ifdef IS_MPI |
682 |
> |
if(worldRank == 0){ |
683 |
> |
#endif |
684 |
> |
//cout << "before resetting the COMVel of sytem is " << zSysCOMVel << endl; |
685 |
> |
#ifdef IS_MPI |
686 |
> |
} |
687 |
> |
#endif |
688 |
> |
|
689 |
> |
// calculate the vz of center of mass of unconstrained molecules and moving z-constrained molecules |
690 |
> |
double MVzOfMovingMols_local; |
691 |
> |
double MVzOfMovingMols; |
692 |
> |
double totalMassOfMovingZMols_local; |
693 |
> |
double totalMassOfMovingZMols; |
694 |
> |
|
695 |
> |
MVzOfMovingMols_local = 0; |
696 |
> |
totalMassOfMovingZMols_local = 0; |
697 |
> |
|
698 |
> |
for(int i =0; i < unconsMols.size(); i++){ |
699 |
> |
unconsMols[i]->getCOMvel(COMvel); |
700 |
> |
MVzOfMovingMols_local += massOfUnconsMols[i] * COMvel[whichDirection]; |
701 |
> |
} |
702 |
> |
|
703 |
> |
for(int i = 0; i < zconsMols.size(); i++){ |
704 |
> |
if (states[i] == zcsMoving){ |
705 |
> |
zconsMols[i]->getCOMvel(COMvel); |
706 |
> |
MVzOfMovingMols_local += massOfZConsMols[i] * COMvel[whichDirection]; |
707 |
> |
totalMassOfMovingZMols_local += massOfZConsMols[i]; |
708 |
> |
} |
709 |
> |
|
710 |
> |
} |
711 |
> |
|
712 |
> |
#ifndef IS_MPI |
713 |
> |
MVzOfMovingMols = MVzOfMovingMols_local; |
714 |
> |
totalMassOfMovingZMols = totalMassOfMovingZMols_local; |
715 |
> |
#else |
716 |
> |
MPI_Allreduce(&MVzOfMovingMols_local, &MVzOfMovingMols, 1, MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
717 |
> |
MPI_Allreduce(&totalMassOfMovingZMols_local, &totalMassOfMovingZMols, 1, MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
718 |
> |
#endif |
719 |
> |
|
720 |
> |
double vzOfMovingMols; |
721 |
> |
vzOfMovingMols = MVzOfMovingMols / (totalMassOfUncons + totalMassOfMovingZMols); |
722 |
> |
|
723 |
> |
//modify the velocites of unconstrained molecules |
724 |
> |
Atom** unconsAtoms; |
725 |
> |
for(int i = 0; i < unconsMols.size(); i++){ |
726 |
> |
|
727 |
> |
unconsAtoms = unconsMols[i]->getMyAtoms(); |
728 |
> |
for(int j = 0; j < unconsMols[i]->getNAtoms();j++){ |
729 |
> |
unconsAtoms[j]->getVel(vel); |
730 |
> |
vel[whichDirection] -= vzOfMovingMols; |
731 |
> |
unconsAtoms[j]->setVel(vel); |
732 |
> |
} |
733 |
> |
|
734 |
> |
} |
735 |
> |
|
736 |
> |
//modify the velocities of moving z-constrained molecuels |
737 |
> |
Atom** movingZAtoms; |
738 |
> |
for(int i = 0; i < zconsMols.size(); i++){ |
739 |
> |
|
740 |
> |
if (states[i] ==zcsMoving){ |
741 |
|
|
742 |
< |
zconsMols[i]->getCOM(COM); |
742 |
> |
movingZAtoms = zconsMols[i]->getMyAtoms(); |
743 |
> |
for(int j = 0; j < zconsMols[i]->getNAtoms(); j++){ |
744 |
> |
movingZAtoms[j]->getVel(vel); |
745 |
> |
vel[whichDirection] -= vzOfMovingMols; |
746 |
> |
movingZAtoms[j]->setVel(vel); |
747 |
> |
} |
748 |
|
|
749 |
< |
deltaZ = zsys + refZ[i] - COM[2]; |
750 |
< |
//update z coordinate |
751 |
< |
zconsAtoms = zconsMols[i]->getMyAtoms(); |
752 |
< |
for(int j =0; j < zconsMols[i]->getNAtoms(); j++){ |
753 |
< |
zconsAtoms[j]->setZ(zconsAtoms[j]->getZ() + deltaZ); |
754 |
< |
} |
749 |
> |
} |
750 |
> |
|
751 |
> |
} |
752 |
> |
|
753 |
> |
|
754 |
> |
zSysCOMVel = calcSysCOMVel(); |
755 |
> |
#ifdef IS_MPI |
756 |
> |
if(worldRank == 0){ |
757 |
> |
#endif |
758 |
> |
//cout << "after resetting the COMVel of moving molecules is " << zSysCOMVel << endl; |
759 |
> |
#ifdef IS_MPI |
760 |
> |
} |
761 |
> |
#endif |
762 |
> |
|
763 |
> |
} |
764 |
> |
|
765 |
> |
/** |
766 |
> |
* |
767 |
> |
*/ |
768 |
> |
|
769 |
> |
template<typename T> void ZConstraint<T>::doZconstraintForce(){ |
770 |
> |
|
771 |
> |
Atom** zconsAtoms; |
772 |
> |
double totalFZ; |
773 |
> |
double totalFZ_local; |
774 |
> |
double COMvel[3]; |
775 |
> |
double COM[3]; |
776 |
> |
double force[3]; |
777 |
> |
|
778 |
> |
//constrain the molecules which do not reach the specified positions |
779 |
|
|
780 |
< |
//calculate z constrain force |
781 |
< |
fz[i] = massOfZConsMols[i]* deltaZ / dt2; |
780 |
> |
//Zero Out the force of z-contrained molecules |
781 |
> |
totalFZ_local = 0; |
782 |
> |
|
783 |
> |
//calculate the total z-contrained force of fixed z-contrained molecules |
784 |
> |
|
785 |
> |
//cout << "before zero out z-constraint force on fixed z-constraint molecuels " |
786 |
> |
// << "total force is " << calcTotalForce() << endl; |
787 |
> |
|
788 |
> |
for(int i = 0; i < zconsMols.size(); i++){ |
789 |
|
|
790 |
+ |
if (states[i] == zcsFixed){ |
791 |
+ |
|
792 |
+ |
zconsMols[i]->getCOM(COM); |
793 |
+ |
zconsAtoms = zconsMols[i]->getMyAtoms(); |
794 |
+ |
|
795 |
+ |
fz[i] = 0; |
796 |
+ |
for(int j =0; j < zconsMols[i]->getNAtoms(); j++) { |
797 |
+ |
zconsAtoms[j]->getFrc(force); |
798 |
+ |
fz[i] += force[whichDirection]; |
799 |
+ |
} |
800 |
+ |
totalFZ_local += fz[i]; |
801 |
+ |
|
802 |
+ |
//cout << "Fixed Molecule\tindex: " << indexOfZConsMols[i] |
803 |
+ |
// <<"\tcurrent zpos: " << COM[whichDirection] |
804 |
+ |
// << "\tcurrent fz: " <<fz[i] << endl; |
805 |
+ |
|
806 |
+ |
|
807 |
+ |
} |
808 |
+ |
|
809 |
+ |
} |
810 |
+ |
|
811 |
+ |
//calculate total z-constraint force |
812 |
+ |
#ifdef IS_MPI |
813 |
+ |
MPI_Allreduce(&totalFZ_local, &totalFZ, 1, MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
814 |
+ |
#else |
815 |
+ |
totalFZ = totalFZ_local; |
816 |
+ |
#endif |
817 |
+ |
|
818 |
+ |
|
819 |
+ |
// apply negative to fixed z-constrained molecues; |
820 |
+ |
force[0]= 0; |
821 |
+ |
force[1]= 0; |
822 |
+ |
force[2]= 0; |
823 |
+ |
|
824 |
+ |
for(int i = 0; i < zconsMols.size(); i++){ |
825 |
+ |
|
826 |
+ |
if (states[i] == zcsFixed){ |
827 |
+ |
|
828 |
+ |
int nAtomOfCurZConsMol = zconsMols[i]->getNAtoms(); |
829 |
+ |
zconsAtoms = zconsMols[i]->getMyAtoms(); |
830 |
+ |
|
831 |
+ |
for(int j =0; j < nAtomOfCurZConsMol; j++) { |
832 |
+ |
//force[whichDirection] = -fz[i]/ nAtomOfCurZConsMol; |
833 |
+ |
force[whichDirection] = - forcePolicy->getZFOfFixedZMols(zconsMols[i], zconsAtoms[j], fz[i]); |
834 |
+ |
zconsAtoms[j]->addFrc(force); |
835 |
+ |
} |
836 |
+ |
|
837 |
+ |
} |
838 |
+ |
|
839 |
+ |
} |
840 |
+ |
|
841 |
+ |
//cout << "after zero out z-constraint force on fixed z-constraint molecuels " |
842 |
+ |
// << "total force is " << calcTotalForce() << endl; |
843 |
+ |
|
844 |
+ |
|
845 |
+ |
force[0]= 0; |
846 |
+ |
force[1]= 0; |
847 |
+ |
force[2]= 0; |
848 |
+ |
|
849 |
+ |
//modify the forces of unconstrained molecules |
850 |
+ |
for(int i = 0; i < unconsMols.size(); i++){ |
851 |
+ |
|
852 |
+ |
Atom** unconsAtoms = unconsMols[i]->getMyAtoms(); |
853 |
+ |
|
854 |
+ |
for(int j = 0; j < unconsMols[i]->getNAtoms(); j++){ |
855 |
+ |
//force[whichDirection] = totalFZ / (totNumOfUnconsAtoms + nMovingZAtoms); |
856 |
+ |
force[whichDirection] = forcePolicy->getZFOfMovingMols(unconsAtoms[j],totalFZ); |
857 |
+ |
unconsAtoms[j]->addFrc(force); |
858 |
+ |
} |
859 |
+ |
|
860 |
+ |
} |
861 |
+ |
|
862 |
+ |
//modify the forces of moving z-constrained molecules |
863 |
+ |
for(int i = 0; i < zconsMols.size(); i++) { |
864 |
+ |
if (states[i] == zcsMoving){ |
865 |
+ |
|
866 |
+ |
Atom** movingZAtoms = zconsMols[i]->getMyAtoms(); |
867 |
+ |
|
868 |
+ |
for(int j = 0; j < zconsMols[i]->getNAtoms(); j++){ |
869 |
+ |
//force[whichDirection] = totalFZ / (totNumOfUnconsAtoms + nMovingZAtoms); |
870 |
+ |
force[whichDirection] = forcePolicy->getZFOfMovingMols(movingZAtoms[j],totalFZ); |
871 |
+ |
movingZAtoms[j]->addFrc(force); |
872 |
+ |
} |
873 |
+ |
} |
874 |
|
} |
875 |
+ |
// cout << "after substracting z-constraint force from moving molecuels " |
876 |
+ |
// << "total force is " << calcTotalForce() << endl; |
877 |
|
|
878 |
+ |
} |
879 |
+ |
|
880 |
+ |
/** |
881 |
+ |
* |
882 |
+ |
* |
883 |
+ |
*/ |
884 |
+ |
|
885 |
+ |
template<typename T> void ZConstraint<T>::doHarmonic(){ |
886 |
+ |
double force[3]; |
887 |
+ |
double harmonicU; |
888 |
+ |
double harmonicF; |
889 |
+ |
double COM[3]; |
890 |
+ |
double diff; |
891 |
+ |
double totalFZ_local; |
892 |
+ |
double totalFZ; |
893 |
+ |
|
894 |
+ |
force[0] = 0; |
895 |
+ |
force[1] = 0; |
896 |
+ |
force[2] = 0; |
897 |
+ |
|
898 |
+ |
totalFZ_local = 0; |
899 |
+ |
|
900 |
+ |
for(int i = 0; i < zconsMols.size(); i++) { |
901 |
+ |
|
902 |
+ |
if (states[i] == zcsMoving){ |
903 |
+ |
zconsMols[i]->getCOM(COM); |
904 |
+ |
// cout << "Moving Molecule\tindex: " << indexOfZConsMols[i] |
905 |
+ |
// << "\tcurrent zpos: " << COM[whichDirection] << endl; |
906 |
+ |
|
907 |
+ |
diff = COM[whichDirection] -zPos[i]; |
908 |
+ |
|
909 |
+ |
harmonicU = 0.5 * kz[i] * diff * diff; |
910 |
+ |
info->lrPot += harmonicU; |
911 |
+ |
|
912 |
+ |
harmonicF = - kz[i] * diff; |
913 |
+ |
totalFZ_local += harmonicF; |
914 |
+ |
|
915 |
+ |
//adjust force |
916 |
+ |
|
917 |
+ |
Atom** movingZAtoms = zconsMols[i]->getMyAtoms(); |
918 |
+ |
|
919 |
+ |
for(int j = 0; j < zconsMols[i]->getNAtoms(); j++){ |
920 |
+ |
//force[whichDirection] = harmonicF / zconsMols[i]->getNAtoms(); |
921 |
+ |
force[whichDirection] = forcePolicy->getHFOfFixedZMols(zconsMols[i], movingZAtoms[j], harmonicF); |
922 |
+ |
movingZAtoms[j]->addFrc(force); |
923 |
+ |
} |
924 |
+ |
} |
925 |
+ |
|
926 |
+ |
} |
927 |
+ |
|
928 |
+ |
#ifndef IS_MPI |
929 |
+ |
totalFZ = totalFZ_local; |
930 |
+ |
#else |
931 |
+ |
MPI_Allreduce(&totalFZ_local, &totalFZ, 1, MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
932 |
+ |
#endif |
933 |
+ |
|
934 |
+ |
//cout << "before substracting harmonic force from moving molecuels " |
935 |
+ |
// << "total force is " << calcTotalForce() << endl; |
936 |
+ |
|
937 |
+ |
force[0]= 0; |
938 |
+ |
force[1]= 0; |
939 |
+ |
force[2]= 0; |
940 |
+ |
|
941 |
+ |
//modify the forces of unconstrained molecules |
942 |
+ |
for(int i = 0; i < unconsMols.size(); i++){ |
943 |
+ |
|
944 |
+ |
Atom** unconsAtoms = unconsMols[i]->getMyAtoms(); |
945 |
+ |
|
946 |
+ |
for(int j = 0; j < unconsMols[i]->getNAtoms(); j++){ |
947 |
+ |
//force[whichDirection] = - totalFZ /totNumOfUnconsAtoms; |
948 |
+ |
force[whichDirection] = - forcePolicy->getHFOfUnconsMols(unconsAtoms[j], totalFZ); |
949 |
+ |
unconsAtoms[j]->addFrc(force); |
950 |
+ |
} |
951 |
+ |
} |
952 |
+ |
|
953 |
+ |
//cout << "after substracting harmonic force from moving molecuels " |
954 |
+ |
// << "total force is " << calcTotalForce() << endl; |
955 |
+ |
|
956 |
+ |
} |
957 |
+ |
|
958 |
+ |
/** |
959 |
+ |
* |
960 |
+ |
*/ |
961 |
+ |
|
962 |
+ |
template<typename T> bool ZConstraint<T>::checkZConsState(){ |
963 |
+ |
double COM[3]; |
964 |
+ |
double diff; |
965 |
+ |
|
966 |
+ |
int changed_local; |
967 |
+ |
int changed; |
968 |
+ |
|
969 |
+ |
changed_local = 0; |
970 |
+ |
|
971 |
+ |
for(int i =0; i < zconsMols.size(); i++){ |
972 |
+ |
|
973 |
+ |
zconsMols[i]->getCOM(COM); |
974 |
+ |
diff = fabs(COM[whichDirection] - zPos[i]); |
975 |
+ |
if ( diff <= zconsTol && states[i] == zcsMoving){ |
976 |
+ |
states[i] = zcsFixed; |
977 |
+ |
changed_local = 1; |
978 |
+ |
} |
979 |
+ |
else if ( diff > zconsTol && states[i] == zcsFixed){ |
980 |
+ |
states[i] = zcsMoving; |
981 |
+ |
changed_local = 1; |
982 |
+ |
} |
983 |
+ |
|
984 |
+ |
} |
985 |
+ |
|
986 |
+ |
#ifndef IS_MPI |
987 |
+ |
changed =changed_local; |
988 |
+ |
#else |
989 |
+ |
MPI_Allreduce(&changed_local, &changed, 1, MPI_INT,MPI_SUM, MPI_COMM_WORLD); |
990 |
+ |
#endif |
991 |
+ |
|
992 |
+ |
return (changed > 0); |
993 |
+ |
|
994 |
+ |
} |
995 |
+ |
|
996 |
+ |
template<typename T> bool ZConstraint<T>::haveFixedZMols(){ |
997 |
+ |
|
998 |
+ |
int havingFixed_local; |
999 |
+ |
int havingFixed; |
1000 |
+ |
|
1001 |
+ |
havingFixed_local = 0; |
1002 |
+ |
|
1003 |
+ |
for(int i = 0; i < zconsMols.size(); i++) |
1004 |
+ |
if (states[i] == zcsFixed){ |
1005 |
+ |
havingFixed_local = 1; |
1006 |
+ |
break; |
1007 |
+ |
} |
1008 |
+ |
|
1009 |
+ |
#ifndef IS_MPI |
1010 |
+ |
havingFixed = havingFixed_local; |
1011 |
+ |
#else |
1012 |
+ |
MPI_Allreduce(&havingFixed_local, &havingFixed, 1, MPI_INT,MPI_SUM, MPI_COMM_WORLD); |
1013 |
+ |
#endif |
1014 |
+ |
|
1015 |
+ |
return (havingFixed > 0); |
1016 |
+ |
} |
1017 |
+ |
|
1018 |
+ |
|
1019 |
+ |
/** |
1020 |
+ |
* |
1021 |
+ |
*/ |
1022 |
+ |
template<typename T> bool ZConstraint<T>::haveMovingZMols(){ |
1023 |
+ |
|
1024 |
+ |
int havingMoving_local; |
1025 |
+ |
int havingMoving; |
1026 |
+ |
|
1027 |
+ |
havingMoving_local = 0; |
1028 |
+ |
|
1029 |
+ |
for(int i = 0; i < zconsMols.size(); i++) |
1030 |
+ |
if (states[i] == zcsMoving){ |
1031 |
+ |
havingMoving_local = 1; |
1032 |
+ |
break; |
1033 |
+ |
} |
1034 |
+ |
|
1035 |
+ |
#ifndef IS_MPI |
1036 |
+ |
havingMoving = havingMoving_local; |
1037 |
+ |
#else |
1038 |
+ |
MPI_Allreduce(&havingMoving_local, &havingMoving, 1, MPI_INT,MPI_SUM, MPI_COMM_WORLD); |
1039 |
+ |
#endif |
1040 |
+ |
|
1041 |
+ |
return (havingMoving > 0); |
1042 |
+ |
|
1043 |
+ |
} |
1044 |
+ |
|
1045 |
+ |
/** |
1046 |
+ |
* |
1047 |
+ |
*/ |
1048 |
+ |
|
1049 |
+ |
template<typename T> double ZConstraint<T>::calcMovingMolsCOMVel() |
1050 |
+ |
{ |
1051 |
+ |
double MVzOfMovingMols_local; |
1052 |
+ |
double MVzOfMovingMols; |
1053 |
+ |
double totalMassOfMovingZMols_local; |
1054 |
+ |
double totalMassOfMovingZMols; |
1055 |
+ |
double COMvel[3]; |
1056 |
|
|
1057 |
+ |
MVzOfMovingMols_local = 0; |
1058 |
+ |
totalMassOfMovingZMols_local = 0; |
1059 |
+ |
|
1060 |
+ |
for(int i =0; i < unconsMols.size(); i++){ |
1061 |
+ |
unconsMols[i]->getCOMvel(COMvel); |
1062 |
+ |
MVzOfMovingMols_local += massOfUnconsMols[i] * COMvel[whichDirection]; |
1063 |
+ |
} |
1064 |
+ |
|
1065 |
+ |
for(int i = 0; i < zconsMols.size(); i++){ |
1066 |
+ |
|
1067 |
+ |
if (states[i] == zcsMoving){ |
1068 |
+ |
zconsMols[i]->getCOMvel(COMvel); |
1069 |
+ |
MVzOfMovingMols_local += massOfZConsMols[i] * COMvel[whichDirection]; |
1070 |
+ |
totalMassOfMovingZMols_local += massOfZConsMols[i]; |
1071 |
+ |
} |
1072 |
+ |
|
1073 |
+ |
} |
1074 |
+ |
|
1075 |
+ |
#ifndef IS_MPI |
1076 |
+ |
MVzOfMovingMols = MVzOfMovingMols_local; |
1077 |
+ |
totalMassOfMovingZMols = totalMassOfMovingZMols_local; |
1078 |
+ |
#else |
1079 |
+ |
MPI_Allreduce(&MVzOfMovingMols_local, &MVzOfMovingMols, 1, MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
1080 |
+ |
MPI_Allreduce(&totalMassOfMovingZMols_local, &totalMassOfMovingZMols, 1, MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
1081 |
+ |
#endif |
1082 |
+ |
|
1083 |
+ |
double vzOfMovingMols; |
1084 |
+ |
vzOfMovingMols = MVzOfMovingMols / (totalMassOfUncons + totalMassOfMovingZMols); |
1085 |
+ |
|
1086 |
+ |
return vzOfMovingMols; |
1087 |
|
} |
1088 |
+ |
|
1089 |
+ |
/** |
1090 |
+ |
* |
1091 |
+ |
*/ |
1092 |
+ |
|
1093 |
+ |
template<typename T> double ZConstraint<T>::calcSysCOMVel() |
1094 |
+ |
{ |
1095 |
+ |
double COMvel[3]; |
1096 |
+ |
double tempMVz_local; |
1097 |
+ |
double tempMVz; |
1098 |
+ |
double massOfZCons_local; |
1099 |
+ |
double massOfZCons; |
1100 |
+ |
|
1101 |
+ |
|
1102 |
+ |
tempMVz_local = 0; |
1103 |
+ |
|
1104 |
+ |
for(int i =0 ; i < nMols; i++){ |
1105 |
+ |
molecules[i].getCOMvel(COMvel); |
1106 |
+ |
tempMVz_local += molecules[i].getTotalMass()*COMvel[whichDirection]; |
1107 |
+ |
} |
1108 |
+ |
|
1109 |
+ |
massOfZCons_local = 0; |
1110 |
+ |
|
1111 |
+ |
for(int i = 0; i < massOfZConsMols.size(); i++){ |
1112 |
+ |
massOfZCons_local += massOfZConsMols[i]; |
1113 |
+ |
} |
1114 |
+ |
#ifndef IS_MPI |
1115 |
+ |
massOfZCons = massOfZCons_local; |
1116 |
+ |
tempMVz = tempMVz_local; |
1117 |
+ |
#else |
1118 |
+ |
MPI_Allreduce(&massOfZCons_local, &massOfZCons, 1, MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
1119 |
+ |
MPI_Allreduce(&tempMVz_local, &tempMVz, 1, MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
1120 |
+ |
#endif |
1121 |
+ |
|
1122 |
+ |
return tempMVz /(totalMassOfUncons + massOfZCons); |
1123 |
+ |
} |
1124 |
+ |
|
1125 |
+ |
/** |
1126 |
+ |
* |
1127 |
+ |
*/ |
1128 |
+ |
|
1129 |
+ |
template<typename T> double ZConstraint<T>::calcTotalForce(){ |
1130 |
+ |
|
1131 |
+ |
double force[3]; |
1132 |
+ |
double totalForce_local; |
1133 |
+ |
double totalForce; |
1134 |
+ |
|
1135 |
+ |
totalForce_local = 0; |
1136 |
+ |
|
1137 |
+ |
for(int i = 0; i < nAtoms; i++){ |
1138 |
+ |
atoms[i]->getFrc(force); |
1139 |
+ |
totalForce_local += force[whichDirection]; |
1140 |
+ |
} |
1141 |
+ |
|
1142 |
+ |
#ifndef IS_MPI |
1143 |
+ |
totalForce = totalForce_local; |
1144 |
+ |
#else |
1145 |
+ |
MPI_Allreduce(&totalForce_local, &totalForce, 1, MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
1146 |
+ |
#endif |
1147 |
+ |
|
1148 |
+ |
return totalForce; |
1149 |
+ |
|
1150 |
+ |
} |
1151 |
+ |
|
1152 |
+ |
/** |
1153 |
+ |
* |
1154 |
+ |
*/ |
1155 |
+ |
|
1156 |
+ |
template<typename T> void ZConstraint<T>::PolicyByNumber::update(){ |
1157 |
+ |
//calculate the number of atoms of moving z-constrained molecules |
1158 |
+ |
int nMovingZAtoms_local; |
1159 |
+ |
int nMovingZAtoms; |
1160 |
+ |
|
1161 |
+ |
nMovingZAtoms_local = 0; |
1162 |
+ |
for(int i = 0; i < (zconsIntegrator->zconsMols).size(); i++) |
1163 |
+ |
if((zconsIntegrator->states)[i] == (zconsIntegrator->zcsMoving)) |
1164 |
+ |
nMovingZAtoms_local += (zconsIntegrator->zconsMols)[i]->getNAtoms(); |
1165 |
+ |
|
1166 |
+ |
#ifdef IS_MPI |
1167 |
+ |
MPI_Allreduce(&nMovingZAtoms_local, &nMovingZAtoms, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD); |
1168 |
+ |
#else |
1169 |
+ |
nMovingZAtoms = nMovingZAtoms_local; |
1170 |
+ |
#endif |
1171 |
+ |
totNumOfMovingAtoms = nMovingZAtoms + zconsIntegrator->totNumOfUnconsAtoms; |
1172 |
+ |
} |
1173 |
+ |
|
1174 |
+ |
template<typename T>double ZConstraint<T>::PolicyByNumber::getZFOfFixedZMols(Molecule* mol, Atom* atom, double totalForce){ |
1175 |
+ |
return totalForce / mol->getNAtoms(); |
1176 |
+ |
} |
1177 |
+ |
|
1178 |
+ |
template<typename T> double ZConstraint<T>::PolicyByNumber::getZFOfMovingMols(Atom* atom, double totalForce){ |
1179 |
+ |
return totalForce / totNumOfMovingAtoms; |
1180 |
+ |
} |
1181 |
+ |
|
1182 |
+ |
template<typename T> double ZConstraint<T>::PolicyByNumber::getHFOfFixedZMols(Molecule* mol, Atom* atom, double totalForce){ |
1183 |
+ |
return totalForce / mol->getNAtoms(); |
1184 |
+ |
} |
1185 |
+ |
|
1186 |
+ |
template<typename T> double ZConstraint<T>::PolicyByNumber::getHFOfUnconsMols(Atom* atom, double totalForce){ |
1187 |
+ |
return totalForce / zconsIntegrator->totNumOfUnconsAtoms; |
1188 |
+ |
} |
1189 |
+ |
|
1190 |
+ |
/** |
1191 |
+ |
* |
1192 |
+ |
*/ |
1193 |
+ |
|
1194 |
+ |
template<typename T> void ZConstraint<T>::PolicyByMass::update(){ |
1195 |
+ |
//calculate the number of atoms of moving z-constrained molecules |
1196 |
+ |
double massOfMovingZAtoms_local; |
1197 |
+ |
double massOfMovingZAtoms; |
1198 |
+ |
|
1199 |
+ |
massOfMovingZAtoms_local = 0; |
1200 |
+ |
for(int i = 0; i < (zconsIntegrator->zconsMols).size(); i++) |
1201 |
+ |
if((zconsIntegrator->states)[i] == (zconsIntegrator->zcsMoving)) |
1202 |
+ |
massOfMovingZAtoms_local += (zconsIntegrator->zconsMols)[i]->getTotalMass(); |
1203 |
+ |
|
1204 |
+ |
#ifdef IS_MPI |
1205 |
+ |
MPI_Allreduce(&massOfMovingZAtoms_local, &massOfMovingZAtoms, 1, MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD); |
1206 |
+ |
#else |
1207 |
+ |
massOfMovingZAtoms = massOfMovingZAtoms_local; |
1208 |
+ |
#endif |
1209 |
+ |
totMassOfMovingAtoms = massOfMovingZAtoms + zconsIntegrator->totalMassOfUncons; |
1210 |
+ |
} |
1211 |
+ |
|
1212 |
+ |
template<typename T> double ZConstraint<T>::PolicyByMass::getZFOfFixedZMols(Molecule* mol, Atom* atom, double totalForce){ |
1213 |
+ |
return totalForce * atom->getMass() / mol->getTotalMass(); |
1214 |
+ |
} |
1215 |
+ |
|
1216 |
+ |
template<typename T> double ZConstraint<T>::PolicyByMass::getZFOfMovingMols( Atom* atom, double totalForce){ |
1217 |
+ |
return totalForce * atom->getMass() / totMassOfMovingAtoms; |
1218 |
+ |
} |
1219 |
+ |
|
1220 |
+ |
template<typename T> double ZConstraint<T>::PolicyByMass::getHFOfFixedZMols(Molecule* mol, Atom* atom, double totalForce){ |
1221 |
+ |
return totalForce * atom->getMass() / mol->getTotalMass(); |
1222 |
+ |
} |
1223 |
+ |
|
1224 |
+ |
template<typename T> double ZConstraint<T>::PolicyByMass::getHFOfUnconsMols(Atom* atom, double totalForce){ |
1225 |
+ |
return totalForce * atom->getMass() / zconsIntegrator->totalMassOfUncons; |
1226 |
+ |
} |
1227 |
+ |
|