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