1 |
gezelter |
507 |
/* |
2 |
gezelter |
246 |
* Copyright (c) 2005 The University of Notre Dame. All Rights Reserved. |
3 |
|
|
* |
4 |
|
|
* The University of Notre Dame grants you ("Licensee") a |
5 |
|
|
* non-exclusive, royalty free, license to use, modify and |
6 |
|
|
* redistribute this software in source and binary code form, provided |
7 |
|
|
* that the following conditions are met: |
8 |
|
|
* |
9 |
gezelter |
1390 |
* 1. Redistributions of source code must retain the above copyright |
10 |
gezelter |
246 |
* notice, this list of conditions and the following disclaimer. |
11 |
|
|
* |
12 |
gezelter |
1390 |
* 2. Redistributions in binary form must reproduce the above copyright |
13 |
gezelter |
246 |
* notice, this list of conditions and the following disclaimer in the |
14 |
|
|
* documentation and/or other materials provided with the |
15 |
|
|
* distribution. |
16 |
|
|
* |
17 |
|
|
* This software is provided "AS IS," without a warranty of any |
18 |
|
|
* kind. All express or implied conditions, representations and |
19 |
|
|
* warranties, including any implied warranty of merchantability, |
20 |
|
|
* fitness for a particular purpose or non-infringement, are hereby |
21 |
|
|
* excluded. The University of Notre Dame and its licensors shall not |
22 |
|
|
* be liable for any damages suffered by licensee as a result of |
23 |
|
|
* using, modifying or distributing the software or its |
24 |
|
|
* derivatives. In no event will the University of Notre Dame or its |
25 |
|
|
* licensors be liable for any lost revenue, profit or data, or for |
26 |
|
|
* direct, indirect, special, consequential, incidental or punitive |
27 |
|
|
* damages, however caused and regardless of the theory of liability, |
28 |
|
|
* arising out of the use of or inability to use software, even if the |
29 |
|
|
* University of Notre Dame has been advised of the possibility of |
30 |
|
|
* such damages. |
31 |
gezelter |
1390 |
* |
32 |
|
|
* SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your |
33 |
|
|
* research, please cite the appropriate papers when you publish your |
34 |
|
|
* work. Good starting points are: |
35 |
|
|
* |
36 |
|
|
* [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005). |
37 |
|
|
* [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006). |
38 |
|
|
* [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008). |
39 |
|
|
* [4] Vardeman & Gezelter, in progress (2009). |
40 |
gezelter |
246 |
*/ |
41 |
|
|
|
42 |
|
|
#include <cmath> |
43 |
|
|
#include "constraints/ZconstraintForceManager.hpp" |
44 |
|
|
#include "integrators/Integrator.hpp" |
45 |
|
|
#include "utils/simError.h" |
46 |
gezelter |
1390 |
#include "utils/PhysicalConstants.hpp" |
47 |
gezelter |
246 |
#include "utils/StringUtils.hpp" |
48 |
gezelter |
1627 |
#ifdef IS_MPI |
49 |
|
|
#include <mpi.h> |
50 |
|
|
#endif |
51 |
|
|
|
52 |
gezelter |
1390 |
namespace OpenMD { |
53 |
gezelter |
507 |
ZconstraintForceManager::ZconstraintForceManager(SimInfo* info): ForceManager(info), infiniteTime(1e31) { |
54 |
gezelter |
246 |
currSnapshot_ = info_->getSnapshotManager()->getCurrentSnapshot(); |
55 |
|
|
Globals* simParam = info_->getSimParams(); |
56 |
|
|
|
57 |
|
|
if (simParam->haveDt()){ |
58 |
gezelter |
507 |
dt_ = simParam->getDt(); |
59 |
gezelter |
246 |
} else { |
60 |
gezelter |
507 |
sprintf(painCave.errMsg, |
61 |
|
|
"Integrator Error: dt is not set\n"); |
62 |
|
|
painCave.isFatal = 1; |
63 |
|
|
simError(); |
64 |
gezelter |
246 |
} |
65 |
|
|
|
66 |
tim |
665 |
if (simParam->haveZconsTime()){ |
67 |
gezelter |
507 |
zconsTime_ = simParam->getZconsTime(); |
68 |
gezelter |
246 |
} |
69 |
|
|
else{ |
70 |
gezelter |
507 |
sprintf(painCave.errMsg, |
71 |
|
|
"ZConstraint error: If you use a ZConstraint,\n" |
72 |
|
|
"\tyou must set zconsTime.\n"); |
73 |
|
|
painCave.isFatal = 1; |
74 |
|
|
simError(); |
75 |
gezelter |
246 |
} |
76 |
|
|
|
77 |
|
|
if (simParam->haveZconsTol()){ |
78 |
gezelter |
507 |
zconsTol_ = simParam->getZconsTol(); |
79 |
gezelter |
246 |
} |
80 |
|
|
else{ |
81 |
gezelter |
507 |
zconsTol_ = 0.01; |
82 |
|
|
sprintf(painCave.errMsg, |
83 |
|
|
"ZConstraint Warning: Tolerance for z-constraint method is not specified.\n" |
84 |
gezelter |
1390 |
"\tOpenMD will use a default value of %f.\n" |
85 |
gezelter |
507 |
"\tTo set the tolerance, use the zconsTol variable.\n", |
86 |
|
|
zconsTol_); |
87 |
|
|
painCave.isFatal = 0; |
88 |
|
|
simError(); |
89 |
gezelter |
246 |
} |
90 |
|
|
|
91 |
|
|
//set zcons gap |
92 |
tim |
665 |
if (simParam->haveZconsGap()){ |
93 |
gezelter |
507 |
usingZconsGap_ = true; |
94 |
|
|
zconsGap_ = simParam->getZconsGap(); |
95 |
gezelter |
246 |
}else { |
96 |
gezelter |
507 |
usingZconsGap_ = false; |
97 |
|
|
zconsGap_ = 0.0; |
98 |
gezelter |
246 |
} |
99 |
|
|
|
100 |
|
|
//set zcons fixtime |
101 |
tim |
665 |
if (simParam->haveZconsFixtime()){ |
102 |
gezelter |
507 |
zconsFixingTime_ = simParam->getZconsFixtime(); |
103 |
gezelter |
246 |
} else { |
104 |
gezelter |
507 |
zconsFixingTime_ = infiniteTime; |
105 |
gezelter |
246 |
} |
106 |
|
|
|
107 |
|
|
//set zconsUsingSMD |
108 |
tim |
665 |
if (simParam->haveZconsUsingSMD()){ |
109 |
gezelter |
507 |
usingSMD_ = simParam->getZconsUsingSMD(); |
110 |
gezelter |
246 |
}else { |
111 |
gezelter |
507 |
usingSMD_ =false; |
112 |
gezelter |
246 |
} |
113 |
|
|
|
114 |
|
|
zconsOutput_ = getPrefix(info_->getFinalConfigFileName()) + ".fz"; |
115 |
|
|
|
116 |
|
|
//estimate the force constant of harmonical potential |
117 |
|
|
Mat3x3d hmat = currSnapshot_->getHmat(); |
118 |
tim |
963 |
RealType halfOfLargestBox = std::max(hmat(0, 0), std::max(hmat(1, 1), hmat(2, 2))) /2; |
119 |
|
|
RealType targetTemp; |
120 |
gezelter |
246 |
if (simParam->haveTargetTemp()) { |
121 |
gezelter |
507 |
targetTemp = simParam->getTargetTemp(); |
122 |
gezelter |
246 |
} else { |
123 |
gezelter |
507 |
targetTemp = 298.0; |
124 |
gezelter |
246 |
} |
125 |
gezelter |
1390 |
RealType zforceConstant = PhysicalConstants::kb * targetTemp / (halfOfLargestBox * halfOfLargestBox); |
126 |
gezelter |
246 |
|
127 |
tim |
770 |
int nZconstraints = simParam->getNZconsStamps(); |
128 |
|
|
std::vector<ZConsStamp*> stamp = simParam->getZconsStamps(); |
129 |
gezelter |
246 |
// |
130 |
|
|
for (int i = 0; i < nZconstraints; i++){ |
131 |
|
|
|
132 |
gezelter |
507 |
ZconstraintParam param; |
133 |
|
|
int zmolIndex = stamp[i]->getMolIndex(); |
134 |
|
|
if (stamp[i]->haveZpos()) { |
135 |
|
|
param.zTargetPos = stamp[i]->getZpos(); |
136 |
|
|
} else { |
137 |
|
|
param.zTargetPos = getZTargetPos(zmolIndex); |
138 |
|
|
} |
139 |
gezelter |
246 |
|
140 |
gezelter |
507 |
param.kz = zforceConstant * stamp[i]->getKratio(); |
141 |
gezelter |
246 |
|
142 |
gezelter |
507 |
if (stamp[i]->haveCantVel()) { |
143 |
|
|
param.cantVel = stamp[i]->getCantVel(); |
144 |
|
|
} else { |
145 |
|
|
param.cantVel = 0.0; |
146 |
|
|
} |
147 |
gezelter |
246 |
|
148 |
gezelter |
507 |
allZMolIndices_.insert(std::make_pair(zmolIndex, param)); |
149 |
gezelter |
246 |
} |
150 |
|
|
|
151 |
|
|
//create fixedMols_, movingMols_ and unconsMols lists |
152 |
|
|
update(); |
153 |
|
|
|
154 |
|
|
//calculate masss of unconstraint molecules in the whole system (never change during the simulation) |
155 |
tim |
963 |
RealType totMassUnconsMols_local = 0.0; |
156 |
gezelter |
246 |
std::vector<Molecule*>::iterator j; |
157 |
|
|
for ( j = unzconsMols_.begin(); j != unzconsMols_.end(); ++j) { |
158 |
gezelter |
507 |
totMassUnconsMols_local += (*j)->getMass(); |
159 |
gezelter |
246 |
} |
160 |
|
|
#ifndef IS_MPI |
161 |
|
|
totMassUnconsMols_ = totMassUnconsMols_local; |
162 |
|
|
#else |
163 |
tim |
963 |
MPI_Allreduce(&totMassUnconsMols_local, &totMassUnconsMols_, 1, MPI_REALTYPE, |
164 |
gezelter |
507 |
MPI_SUM, MPI_COMM_WORLD); |
165 |
gezelter |
246 |
#endif |
166 |
|
|
|
167 |
|
|
// creat zconsWriter |
168 |
|
|
fzOut = new ZConsWriter(info_, zconsOutput_.c_str()); |
169 |
|
|
|
170 |
|
|
if (!fzOut){ |
171 |
|
|
sprintf(painCave.errMsg, "Fail to create ZConsWriter\n"); |
172 |
|
|
painCave.isFatal = 1; |
173 |
|
|
simError(); |
174 |
|
|
} |
175 |
|
|
|
176 |
gezelter |
507 |
} |
177 |
gezelter |
246 |
|
178 |
gezelter |
507 |
ZconstraintForceManager::~ZconstraintForceManager(){ |
179 |
gezelter |
246 |
|
180 |
gezelter |
507 |
if (fzOut){ |
181 |
|
|
delete fzOut; |
182 |
|
|
} |
183 |
|
|
|
184 |
gezelter |
246 |
} |
185 |
|
|
|
186 |
gezelter |
507 |
void ZconstraintForceManager::update(){ |
187 |
gezelter |
246 |
fixedZMols_.clear(); |
188 |
|
|
movingZMols_.clear(); |
189 |
|
|
unzconsMols_.clear(); |
190 |
|
|
|
191 |
|
|
for (std::map<int, ZconstraintParam>::iterator i = allZMolIndices_.begin(); i != allZMolIndices_.end(); ++i) { |
192 |
|
|
#ifdef IS_MPI |
193 |
gezelter |
507 |
if (info_->getMolToProc(i->first) == worldRank) { |
194 |
gezelter |
246 |
#endif |
195 |
gezelter |
507 |
ZconstraintMol zmol; |
196 |
|
|
zmol.mol = info_->getMoleculeByGlobalIndex(i->first); |
197 |
|
|
assert(zmol.mol); |
198 |
|
|
zmol.param = i->second; |
199 |
|
|
zmol.cantPos = zmol.param.zTargetPos; /**@todo fixed me when zmol migrate, it is incorrect*/ |
200 |
|
|
Vector3d com = zmol.mol->getCom(); |
201 |
tim |
963 |
RealType diff = fabs(zmol.param.zTargetPos - com[whichDirection]); |
202 |
gezelter |
507 |
if (diff < zconsTol_) { |
203 |
|
|
fixedZMols_.push_back(zmol); |
204 |
|
|
} else { |
205 |
|
|
movingZMols_.push_back(zmol); |
206 |
|
|
} |
207 |
gezelter |
246 |
|
208 |
|
|
#ifdef IS_MPI |
209 |
gezelter |
507 |
} |
210 |
gezelter |
246 |
#endif |
211 |
|
|
} |
212 |
|
|
|
213 |
|
|
calcTotalMassMovingZMols(); |
214 |
|
|
|
215 |
|
|
std::set<int> zmolSet; |
216 |
|
|
for (std::list<ZconstraintMol>::iterator i = movingZMols_.begin(); i != movingZMols_.end(); ++i) { |
217 |
gezelter |
507 |
zmolSet.insert(i->mol->getGlobalIndex()); |
218 |
gezelter |
246 |
} |
219 |
|
|
|
220 |
|
|
for (std::list<ZconstraintMol>::iterator i = fixedZMols_.begin(); i != fixedZMols_.end(); ++i) { |
221 |
gezelter |
507 |
zmolSet.insert(i->mol->getGlobalIndex()); |
222 |
gezelter |
246 |
} |
223 |
|
|
|
224 |
|
|
SimInfo::MoleculeIterator mi; |
225 |
|
|
Molecule* mol; |
226 |
|
|
for(mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) { |
227 |
gezelter |
507 |
if (zmolSet.find(mol->getGlobalIndex()) == zmolSet.end()) { |
228 |
|
|
unzconsMols_.push_back(mol); |
229 |
|
|
} |
230 |
gezelter |
246 |
} |
231 |
|
|
|
232 |
gezelter |
507 |
} |
233 |
gezelter |
246 |
|
234 |
gezelter |
507 |
bool ZconstraintForceManager::isZMol(Molecule* mol){ |
235 |
gezelter |
246 |
return allZMolIndices_.find(mol->getGlobalIndex()) == allZMolIndices_.end() ? false : true; |
236 |
gezelter |
507 |
} |
237 |
gezelter |
246 |
|
238 |
gezelter |
507 |
void ZconstraintForceManager::init(){ |
239 |
gezelter |
246 |
|
240 |
gezelter |
507 |
//zero out the velocities of center of mass of unconstrained molecules |
241 |
|
|
//and the velocities of center of mass of every single z-constrained molecueles |
242 |
|
|
zeroVelocity(); |
243 |
gezelter |
246 |
|
244 |
gezelter |
507 |
currZconsTime_ = currSnapshot_->getTime(); |
245 |
|
|
} |
246 |
gezelter |
246 |
|
247 |
gezelter |
1464 |
void ZconstraintForceManager::calcForces(){ |
248 |
|
|
ForceManager::calcForces(); |
249 |
gezelter |
246 |
|
250 |
|
|
if (usingZconsGap_){ |
251 |
gezelter |
507 |
updateZPos(); |
252 |
gezelter |
246 |
} |
253 |
|
|
|
254 |
|
|
if (checkZConsState()){ |
255 |
gezelter |
507 |
zeroVelocity(); |
256 |
|
|
calcTotalMassMovingZMols(); |
257 |
gezelter |
246 |
} |
258 |
|
|
|
259 |
|
|
//do zconstraint force; |
260 |
|
|
if (haveFixedZMols()){ |
261 |
gezelter |
507 |
doZconstraintForce(); |
262 |
gezelter |
246 |
} |
263 |
|
|
|
264 |
|
|
//use external force to move the molecules to the specified positions |
265 |
|
|
if (haveMovingZMols()){ |
266 |
gezelter |
507 |
doHarmonic(); |
267 |
gezelter |
246 |
} |
268 |
|
|
|
269 |
|
|
//write out forces and current positions of z-constraint molecules |
270 |
|
|
if (currSnapshot_->getTime() >= currZconsTime_){ |
271 |
gezelter |
507 |
std::list<ZconstraintMol>::iterator i; |
272 |
|
|
Vector3d com; |
273 |
|
|
for(i = fixedZMols_.begin(); i != fixedZMols_.end(); ++i) { |
274 |
|
|
com = i->mol->getCom(); |
275 |
|
|
i->zpos = com[whichDirection]; |
276 |
|
|
} |
277 |
gezelter |
246 |
|
278 |
gezelter |
507 |
fzOut->writeFZ(fixedZMols_); |
279 |
|
|
currZconsTime_ += zconsTime_; |
280 |
gezelter |
246 |
} |
281 |
gezelter |
507 |
} |
282 |
gezelter |
246 |
|
283 |
gezelter |
507 |
void ZconstraintForceManager::zeroVelocity(){ |
284 |
gezelter |
246 |
|
285 |
|
|
Vector3d comVel; |
286 |
|
|
Vector3d vel; |
287 |
|
|
std::list<ZconstraintMol>::iterator i; |
288 |
|
|
Molecule* mol; |
289 |
|
|
StuntDouble* integrableObject; |
290 |
|
|
Molecule::IntegrableObjectIterator ii; |
291 |
|
|
|
292 |
|
|
//zero out the velocities of center of mass of fixed z-constrained molecules |
293 |
|
|
for(i = fixedZMols_.begin(); i != fixedZMols_.end(); ++i) { |
294 |
gezelter |
507 |
mol = i->mol; |
295 |
|
|
comVel = mol->getComVel(); |
296 |
|
|
for(integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
297 |
|
|
integrableObject = mol->nextIntegrableObject(ii)) { |
298 |
|
|
vel = integrableObject->getVel(); |
299 |
|
|
vel[whichDirection] -= comVel[whichDirection]; |
300 |
|
|
integrableObject->setVel(vel); |
301 |
|
|
} |
302 |
gezelter |
246 |
} |
303 |
|
|
|
304 |
|
|
// calculate the vz of center of mass of moving molecules(include unconstrained molecules |
305 |
|
|
// and moving z-constrained molecules) |
306 |
tim |
963 |
RealType pzMovingMols_local = 0.0; |
307 |
|
|
RealType pzMovingMols; |
308 |
gezelter |
246 |
|
309 |
|
|
for ( i = movingZMols_.begin(); i != movingZMols_.end(); ++i) { |
310 |
gezelter |
507 |
mol = i->mol; |
311 |
|
|
comVel = mol->getComVel(); |
312 |
|
|
pzMovingMols_local += mol->getMass() * comVel[whichDirection]; |
313 |
gezelter |
246 |
} |
314 |
|
|
|
315 |
|
|
std::vector<Molecule*>::iterator j; |
316 |
|
|
for ( j = unzconsMols_.begin(); j != unzconsMols_.end(); ++j) { |
317 |
gezelter |
507 |
mol =*j; |
318 |
|
|
comVel = mol->getComVel(); |
319 |
|
|
pzMovingMols_local += mol->getMass() * comVel[whichDirection]; |
320 |
gezelter |
246 |
} |
321 |
|
|
|
322 |
|
|
#ifndef IS_MPI |
323 |
|
|
pzMovingMols = pzMovingMols_local; |
324 |
|
|
#else |
325 |
tim |
963 |
MPI_Allreduce(&pzMovingMols_local, &pzMovingMols, 1, MPI_REALTYPE, |
326 |
gezelter |
507 |
MPI_SUM, MPI_COMM_WORLD); |
327 |
gezelter |
246 |
#endif |
328 |
|
|
|
329 |
tim |
963 |
RealType vzMovingMols = pzMovingMols / (totMassMovingZMols_ + totMassUnconsMols_); |
330 |
gezelter |
246 |
|
331 |
|
|
//modify the velocities of moving z-constrained molecuels |
332 |
|
|
for ( i = movingZMols_.begin(); i != movingZMols_.end(); ++i) { |
333 |
gezelter |
507 |
mol = i->mol; |
334 |
|
|
for(integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
335 |
|
|
integrableObject = mol->nextIntegrableObject(ii)) { |
336 |
gezelter |
246 |
|
337 |
gezelter |
507 |
vel = integrableObject->getVel(); |
338 |
|
|
vel[whichDirection] -= vzMovingMols; |
339 |
|
|
integrableObject->setVel(vel); |
340 |
|
|
} |
341 |
gezelter |
246 |
} |
342 |
|
|
|
343 |
|
|
//modify the velocites of unconstrained molecules |
344 |
|
|
for ( j = unzconsMols_.begin(); j != unzconsMols_.end(); ++j) { |
345 |
gezelter |
507 |
mol =*j; |
346 |
|
|
for(integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
347 |
|
|
integrableObject = mol->nextIntegrableObject(ii)) { |
348 |
gezelter |
246 |
|
349 |
gezelter |
507 |
vel = integrableObject->getVel(); |
350 |
|
|
vel[whichDirection] -= vzMovingMols; |
351 |
|
|
integrableObject->setVel(vel); |
352 |
|
|
} |
353 |
gezelter |
246 |
} |
354 |
|
|
|
355 |
gezelter |
507 |
} |
356 |
gezelter |
246 |
|
357 |
|
|
|
358 |
gezelter |
507 |
void ZconstraintForceManager::doZconstraintForce(){ |
359 |
tim |
963 |
RealType totalFZ; |
360 |
|
|
RealType totalFZ_local; |
361 |
gezelter |
507 |
Vector3d com; |
362 |
|
|
Vector3d force(0.0); |
363 |
gezelter |
246 |
|
364 |
gezelter |
507 |
//constrain the molecules which do not reach the specified positions |
365 |
gezelter |
246 |
|
366 |
gezelter |
507 |
//Zero Out the force of z-contrained molecules |
367 |
|
|
totalFZ_local = 0; |
368 |
gezelter |
246 |
|
369 |
|
|
|
370 |
|
|
//calculate the total z-contrained force of fixed z-contrained molecules |
371 |
|
|
std::list<ZconstraintMol>::iterator i; |
372 |
|
|
Molecule* mol; |
373 |
|
|
StuntDouble* integrableObject; |
374 |
|
|
Molecule::IntegrableObjectIterator ii; |
375 |
|
|
|
376 |
|
|
for ( i = fixedZMols_.begin(); i != fixedZMols_.end(); ++i) { |
377 |
gezelter |
507 |
mol = i->mol; |
378 |
|
|
i->fz = 0.0; |
379 |
|
|
for(integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
380 |
|
|
integrableObject = mol->nextIntegrableObject(ii)) { |
381 |
gezelter |
246 |
|
382 |
gezelter |
507 |
force = integrableObject->getFrc(); |
383 |
|
|
i->fz += force[whichDirection]; |
384 |
|
|
} |
385 |
|
|
totalFZ_local += i->fz; |
386 |
gezelter |
246 |
} |
387 |
|
|
|
388 |
gezelter |
507 |
//calculate total z-constraint force |
389 |
gezelter |
246 |
#ifdef IS_MPI |
390 |
tim |
963 |
MPI_Allreduce(&totalFZ_local, &totalFZ, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD); |
391 |
gezelter |
246 |
#else |
392 |
gezelter |
507 |
totalFZ = totalFZ_local; |
393 |
gezelter |
246 |
#endif |
394 |
|
|
|
395 |
|
|
|
396 |
gezelter |
507 |
// apply negative to fixed z-constrained molecues; |
397 |
gezelter |
246 |
for ( i = fixedZMols_.begin(); i != fixedZMols_.end(); ++i) { |
398 |
gezelter |
507 |
mol = i->mol; |
399 |
|
|
for(integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
400 |
|
|
integrableObject = mol->nextIntegrableObject(ii)) { |
401 |
gezelter |
246 |
|
402 |
gezelter |
507 |
force[whichDirection] = -getZFOfFixedZMols(mol, integrableObject, i->fz); |
403 |
|
|
integrableObject->addFrc(force); |
404 |
|
|
} |
405 |
gezelter |
246 |
} |
406 |
|
|
|
407 |
gezelter |
507 |
//modify the forces of moving z-constrained molecules |
408 |
gezelter |
246 |
for ( i = movingZMols_.begin(); i != movingZMols_.end(); ++i) { |
409 |
gezelter |
507 |
mol = i->mol; |
410 |
|
|
for(integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
411 |
|
|
integrableObject = mol->nextIntegrableObject(ii)) { |
412 |
gezelter |
246 |
|
413 |
gezelter |
507 |
force[whichDirection] = -getZFOfMovingMols(mol,totalFZ); |
414 |
|
|
integrableObject->addFrc(force); |
415 |
|
|
} |
416 |
gezelter |
246 |
} |
417 |
|
|
|
418 |
gezelter |
507 |
//modify the forces of unconstrained molecules |
419 |
gezelter |
246 |
std::vector<Molecule*>::iterator j; |
420 |
|
|
for ( j = unzconsMols_.begin(); j != unzconsMols_.end(); ++j) { |
421 |
gezelter |
507 |
mol =*j; |
422 |
|
|
for(integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
423 |
|
|
integrableObject = mol->nextIntegrableObject(ii)) { |
424 |
gezelter |
246 |
|
425 |
gezelter |
507 |
force[whichDirection] = -getZFOfMovingMols(mol, totalFZ); |
426 |
|
|
integrableObject->addFrc(force); |
427 |
|
|
} |
428 |
gezelter |
246 |
} |
429 |
|
|
|
430 |
gezelter |
507 |
} |
431 |
gezelter |
246 |
|
432 |
|
|
|
433 |
gezelter |
507 |
void ZconstraintForceManager::doHarmonic(){ |
434 |
tim |
963 |
RealType totalFZ; |
435 |
gezelter |
246 |
Vector3d force(0.0); |
436 |
|
|
Vector3d com; |
437 |
tim |
963 |
RealType totalFZ_local = 0; |
438 |
gezelter |
246 |
std::list<ZconstraintMol>::iterator i; |
439 |
|
|
StuntDouble* integrableObject; |
440 |
|
|
Molecule::IntegrableObjectIterator ii; |
441 |
|
|
Molecule* mol; |
442 |
|
|
for ( i = movingZMols_.begin(); i != movingZMols_.end(); ++i) { |
443 |
gezelter |
507 |
mol = i->mol; |
444 |
|
|
com = mol->getCom(); |
445 |
tim |
963 |
RealType resPos = usingSMD_? i->cantPos : i->param.zTargetPos; |
446 |
|
|
RealType diff = com[whichDirection] - resPos; |
447 |
|
|
RealType harmonicU = 0.5 * i->param.kz * diff * diff; |
448 |
gezelter |
507 |
currSnapshot_->statData[Stats::LONG_RANGE_POTENTIAL] += harmonicU; |
449 |
tim |
963 |
RealType harmonicF = -i->param.kz * diff; |
450 |
gezelter |
507 |
totalFZ_local += harmonicF; |
451 |
gezelter |
246 |
|
452 |
gezelter |
507 |
//adjust force |
453 |
|
|
for(integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
454 |
|
|
integrableObject = mol->nextIntegrableObject(ii)) { |
455 |
gezelter |
246 |
|
456 |
gezelter |
507 |
force[whichDirection] = getHFOfFixedZMols(mol, integrableObject, harmonicF); |
457 |
|
|
integrableObject->addFrc(force); |
458 |
|
|
} |
459 |
gezelter |
246 |
} |
460 |
|
|
|
461 |
|
|
#ifndef IS_MPI |
462 |
gezelter |
507 |
totalFZ = totalFZ_local; |
463 |
gezelter |
246 |
#else |
464 |
tim |
963 |
MPI_Allreduce(&totalFZ_local, &totalFZ, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD); |
465 |
gezelter |
246 |
#endif |
466 |
|
|
|
467 |
|
|
//modify the forces of unconstrained molecules |
468 |
|
|
std::vector<Molecule*>::iterator j; |
469 |
|
|
for ( j = unzconsMols_.begin(); j != unzconsMols_.end(); ++j) { |
470 |
gezelter |
507 |
mol = *j; |
471 |
|
|
for(integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
472 |
|
|
integrableObject = mol->nextIntegrableObject(ii)) { |
473 |
gezelter |
246 |
|
474 |
gezelter |
507 |
force[whichDirection] = getHFOfUnconsMols(mol, totalFZ); |
475 |
|
|
integrableObject->addFrc(force); |
476 |
|
|
} |
477 |
gezelter |
246 |
} |
478 |
|
|
|
479 |
gezelter |
507 |
} |
480 |
gezelter |
246 |
|
481 |
gezelter |
507 |
bool ZconstraintForceManager::checkZConsState(){ |
482 |
gezelter |
246 |
Vector3d com; |
483 |
tim |
963 |
RealType diff; |
484 |
gezelter |
246 |
int changed_local = 0; |
485 |
|
|
|
486 |
|
|
std::list<ZconstraintMol>::iterator i; |
487 |
|
|
std::list<ZconstraintMol>::iterator j; |
488 |
|
|
|
489 |
|
|
std::list<ZconstraintMol> newMovingZMols; |
490 |
|
|
for ( i = fixedZMols_.begin(); i != fixedZMols_.end();) { |
491 |
gezelter |
507 |
com = i->mol->getCom(); |
492 |
|
|
diff = fabs(com[whichDirection] - i->param.zTargetPos); |
493 |
|
|
if (diff > zconsTol_) { |
494 |
|
|
if (usingZconsGap_) { |
495 |
|
|
i->endFixingTime = infiniteTime; |
496 |
|
|
} |
497 |
|
|
j = i++; |
498 |
|
|
newMovingZMols.push_back(*j); |
499 |
|
|
fixedZMols_.erase(j); |
500 |
gezelter |
246 |
|
501 |
gezelter |
507 |
changed_local = 1; |
502 |
|
|
}else { |
503 |
|
|
++i; |
504 |
|
|
} |
505 |
gezelter |
246 |
} |
506 |
|
|
|
507 |
|
|
std::list<ZconstraintMol> newFixedZMols; |
508 |
|
|
for ( i = movingZMols_.begin(); i != movingZMols_.end();) { |
509 |
gezelter |
507 |
com = i->mol->getCom(); |
510 |
|
|
diff = fabs(com[whichDirection] - i->param.zTargetPos); |
511 |
|
|
if (diff <= zconsTol_) { |
512 |
|
|
if (usingZconsGap_) { |
513 |
|
|
i->endFixingTime = currSnapshot_->getTime() + zconsFixingTime_; |
514 |
|
|
} |
515 |
|
|
//this moving zconstraint molecule is about to fixed |
516 |
|
|
//moved this molecule to |
517 |
|
|
j = i++; |
518 |
|
|
newFixedZMols.push_back(*j); |
519 |
|
|
movingZMols_.erase(j); |
520 |
|
|
changed_local = 1; |
521 |
|
|
}else { |
522 |
|
|
++i; |
523 |
|
|
} |
524 |
gezelter |
246 |
} |
525 |
|
|
|
526 |
|
|
//merge the lists |
527 |
|
|
fixedZMols_.insert(fixedZMols_.end(), newFixedZMols.begin(), newFixedZMols.end()); |
528 |
|
|
movingZMols_.insert(movingZMols_.end(), newMovingZMols.begin(), newMovingZMols.end()); |
529 |
|
|
|
530 |
|
|
int changed; |
531 |
|
|
#ifndef IS_MPI |
532 |
|
|
changed = changed_local; |
533 |
|
|
#else |
534 |
|
|
MPI_Allreduce(&changed_local, &changed, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD); |
535 |
|
|
#endif |
536 |
|
|
|
537 |
|
|
return (changed > 0); |
538 |
gezelter |
507 |
} |
539 |
gezelter |
246 |
|
540 |
gezelter |
507 |
bool ZconstraintForceManager::haveFixedZMols(){ |
541 |
|
|
int havingFixed; |
542 |
|
|
int havingFixed_local = fixedZMols_.empty() ? 0 : 1; |
543 |
gezelter |
246 |
|
544 |
|
|
#ifndef IS_MPI |
545 |
gezelter |
507 |
havingFixed = havingFixed_local; |
546 |
gezelter |
246 |
#else |
547 |
gezelter |
507 |
MPI_Allreduce(&havingFixed_local, &havingFixed, 1, MPI_INT, MPI_SUM, |
548 |
|
|
MPI_COMM_WORLD); |
549 |
gezelter |
246 |
#endif |
550 |
|
|
|
551 |
gezelter |
507 |
return havingFixed > 0; |
552 |
|
|
} |
553 |
gezelter |
246 |
|
554 |
|
|
|
555 |
gezelter |
507 |
bool ZconstraintForceManager::haveMovingZMols(){ |
556 |
|
|
int havingMoving_local; |
557 |
|
|
int havingMoving; |
558 |
gezelter |
246 |
|
559 |
gezelter |
507 |
havingMoving_local = movingZMols_.empty()? 0 : 1; |
560 |
gezelter |
246 |
|
561 |
|
|
#ifndef IS_MPI |
562 |
gezelter |
507 |
havingMoving = havingMoving_local; |
563 |
gezelter |
246 |
#else |
564 |
gezelter |
507 |
MPI_Allreduce(&havingMoving_local, &havingMoving, 1, MPI_INT, MPI_SUM, |
565 |
|
|
MPI_COMM_WORLD); |
566 |
gezelter |
246 |
#endif |
567 |
|
|
|
568 |
gezelter |
507 |
return havingMoving > 0; |
569 |
|
|
} |
570 |
gezelter |
246 |
|
571 |
gezelter |
507 |
void ZconstraintForceManager::calcTotalMassMovingZMols(){ |
572 |
gezelter |
246 |
|
573 |
tim |
963 |
RealType totMassMovingZMols_local = 0.0; |
574 |
gezelter |
246 |
std::list<ZconstraintMol>::iterator i; |
575 |
|
|
for ( i = movingZMols_.begin(); i != movingZMols_.end(); ++i) { |
576 |
gezelter |
507 |
totMassMovingZMols_local += i->mol->getMass(); |
577 |
gezelter |
246 |
} |
578 |
|
|
|
579 |
|
|
#ifdef IS_MPI |
580 |
tim |
963 |
MPI_Allreduce(&totMassMovingZMols_local, &totMassMovingZMols_, 1, MPI_REALTYPE, |
581 |
gezelter |
507 |
MPI_SUM, MPI_COMM_WORLD); |
582 |
gezelter |
246 |
#else |
583 |
|
|
totMassMovingZMols_ = totMassMovingZMols_local; |
584 |
|
|
#endif |
585 |
|
|
|
586 |
gezelter |
507 |
} |
587 |
gezelter |
246 |
|
588 |
tim |
963 |
RealType ZconstraintForceManager::getZFOfFixedZMols(Molecule* mol, StuntDouble* sd, RealType totalForce){ |
589 |
gezelter |
507 |
return totalForce * sd->getMass() / mol->getMass(); |
590 |
|
|
} |
591 |
gezelter |
246 |
|
592 |
tim |
963 |
RealType ZconstraintForceManager::getZFOfMovingMols(Molecule* mol, RealType totalForce){ |
593 |
gezelter |
507 |
return totalForce * mol->getMass() / (totMassUnconsMols_ + totMassMovingZMols_); |
594 |
|
|
} |
595 |
gezelter |
246 |
|
596 |
tim |
963 |
RealType ZconstraintForceManager::getHFOfFixedZMols(Molecule* mol, StuntDouble*sd, RealType totalForce){ |
597 |
gezelter |
507 |
return totalForce * sd->getMass() / mol->getMass(); |
598 |
|
|
} |
599 |
gezelter |
246 |
|
600 |
tim |
963 |
RealType ZconstraintForceManager::getHFOfUnconsMols(Molecule* mol, RealType totalForce){ |
601 |
gezelter |
507 |
return totalForce * mol->getMass() / totMassUnconsMols_; |
602 |
|
|
} |
603 |
gezelter |
246 |
|
604 |
gezelter |
507 |
void ZconstraintForceManager::updateZPos(){ |
605 |
tim |
963 |
RealType curTime = currSnapshot_->getTime(); |
606 |
gezelter |
246 |
std::list<ZconstraintMol>::iterator i; |
607 |
|
|
for ( i = fixedZMols_.begin(); i != fixedZMols_.end(); ++i) { |
608 |
gezelter |
507 |
i->param.zTargetPos += zconsGap_; |
609 |
gezelter |
246 |
} |
610 |
gezelter |
507 |
} |
611 |
gezelter |
246 |
|
612 |
gezelter |
507 |
void ZconstraintForceManager::updateCantPos(){ |
613 |
gezelter |
246 |
std::list<ZconstraintMol>::iterator i; |
614 |
|
|
for ( i = movingZMols_.begin(); i != movingZMols_.end(); ++i) { |
615 |
gezelter |
507 |
i->cantPos += i->param.cantVel * dt_; |
616 |
gezelter |
246 |
} |
617 |
gezelter |
507 |
} |
618 |
gezelter |
246 |
|
619 |
tim |
963 |
RealType ZconstraintForceManager::getZTargetPos(int index){ |
620 |
|
|
RealType zTargetPos; |
621 |
gezelter |
246 |
#ifndef IS_MPI |
622 |
|
|
Molecule* mol = info_->getMoleculeByGlobalIndex(index); |
623 |
|
|
assert(mol); |
624 |
|
|
Vector3d com = mol->getCom(); |
625 |
|
|
zTargetPos = com[whichDirection]; |
626 |
|
|
#else |
627 |
|
|
int whicProc = info_->getMolToProc(index); |
628 |
tim |
963 |
MPI_Bcast(&zTargetPos, 1, MPI_REALTYPE, whicProc, MPI_COMM_WORLD); |
629 |
gezelter |
246 |
#endif |
630 |
|
|
return zTargetPos; |
631 |
gezelter |
507 |
} |
632 |
gezelter |
246 |
|
633 |
|
|
} |