ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/constraints/ZconstraintForceManager.cpp
Revision: 1782
Committed: Wed Aug 22 02:28:28 2012 UTC (12 years, 8 months ago) by gezelter
File size: 19130 byte(s)
Log Message:
MERGE OpenMD development branch 1465:1781 into trunk

File Contents

# Content
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* sd;
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
296 mol = i->mol;
297 comVel = mol->getComVel();
298
299 for(sd = mol->beginIntegrableObject(ii); sd != NULL;
300 sd = mol->nextIntegrableObject(ii)) {
301
302 vel = sd->getVel();
303 vel[whichDirection] -= comVel[whichDirection];
304 sd->setVel(vel);
305 }
306 }
307
308 // calculate the vz of center of mass of moving molecules(include unconstrained molecules
309 // and moving z-constrained molecules)
310 RealType pzMovingMols_local = 0.0;
311 RealType pzMovingMols;
312
313 for ( i = movingZMols_.begin(); i != movingZMols_.end(); ++i) {
314 mol = i->mol;
315 comVel = mol->getComVel();
316 pzMovingMols_local += mol->getMass() * comVel[whichDirection];
317 }
318
319 std::vector<Molecule*>::iterator j;
320 for ( j = unzconsMols_.begin(); j != unzconsMols_.end(); ++j) {
321 mol =*j;
322 comVel = mol->getComVel();
323 pzMovingMols_local += mol->getMass() * comVel[whichDirection];
324 }
325
326 #ifndef IS_MPI
327 pzMovingMols = pzMovingMols_local;
328 #else
329 MPI_Allreduce(&pzMovingMols_local, &pzMovingMols, 1, MPI_REALTYPE,
330 MPI_SUM, MPI_COMM_WORLD);
331 #endif
332
333 RealType vzMovingMols = pzMovingMols / (totMassMovingZMols_ + totMassUnconsMols_);
334
335 //modify the velocities of moving z-constrained molecuels
336 for ( i = movingZMols_.begin(); i != movingZMols_.end(); ++i) {
337
338 mol = i->mol;
339
340 for(sd = mol->beginIntegrableObject(ii); sd != NULL;
341 sd = mol->nextIntegrableObject(ii)) {
342
343 vel = sd->getVel();
344 vel[whichDirection] -= vzMovingMols;
345 sd->setVel(vel);
346 }
347 }
348
349 //modify the velocites of unconstrained molecules
350 for ( j = unzconsMols_.begin(); j != unzconsMols_.end(); ++j) {
351
352 mol =*j;
353
354 for(sd = mol->beginIntegrableObject(ii); sd != NULL;
355 sd = mol->nextIntegrableObject(ii)) {
356
357 vel = sd->getVel();
358 vel[whichDirection] -= vzMovingMols;
359 sd->setVel(vel);
360 }
361 }
362
363 }
364
365
366 void ZconstraintForceManager::doZconstraintForce(){
367 RealType totalFZ;
368 RealType totalFZ_local;
369 Vector3d com;
370 Vector3d force(0.0);
371
372 //constrain the molecules which do not reach the specified positions
373
374 //Zero Out the force of z-contrained molecules
375 totalFZ_local = 0;
376
377
378 //calculate the total z-contrained force of fixed z-contrained molecules
379 std::list<ZconstraintMol>::iterator i;
380 Molecule* mol;
381 StuntDouble* sd;
382 Molecule::IntegrableObjectIterator ii;
383
384 for ( i = fixedZMols_.begin(); i != fixedZMols_.end(); ++i) {
385
386 mol = i->mol;
387 i->fz = 0.0;
388
389 for( sd = mol->beginIntegrableObject(ii); sd != NULL;
390 sd = mol->nextIntegrableObject(ii)) {
391
392 force = sd->getFrc();
393 i->fz += force[whichDirection];
394 }
395 totalFZ_local += i->fz;
396 }
397
398 //calculate total z-constraint force
399 #ifdef IS_MPI
400 MPI_Allreduce(&totalFZ_local, &totalFZ, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
401 #else
402 totalFZ = totalFZ_local;
403 #endif
404
405
406 // apply negative to fixed z-constrained molecues;
407 for ( i = fixedZMols_.begin(); i != fixedZMols_.end(); ++i) {
408
409 mol = i->mol;
410
411 for(sd = mol->beginIntegrableObject(ii); sd != NULL;
412 sd = mol->nextIntegrableObject(ii)) {
413
414 force[whichDirection] = -getZFOfFixedZMols(mol, sd, i->fz);
415 sd->addFrc(force);
416 }
417 }
418
419 //modify the forces of moving z-constrained molecules
420 for ( i = movingZMols_.begin(); i != movingZMols_.end(); ++i) {
421
422 mol = i->mol;
423
424 for(sd = mol->beginIntegrableObject(ii); sd != NULL;
425 sd = mol->nextIntegrableObject(ii)) {
426
427 force[whichDirection] = -getZFOfMovingMols(mol,totalFZ);
428 sd->addFrc(force);
429 }
430 }
431
432 //modify the forces of unconstrained molecules
433 std::vector<Molecule*>::iterator j;
434 for ( j = unzconsMols_.begin(); j != unzconsMols_.end(); ++j) {
435
436 mol =*j;
437
438 for(sd = mol->beginIntegrableObject(ii); sd != NULL;
439 sd = mol->nextIntegrableObject(ii)) {
440
441 force[whichDirection] = -getZFOfMovingMols(mol, totalFZ);
442 sd->addFrc(force);
443 }
444 }
445
446 }
447
448
449 void ZconstraintForceManager::doHarmonic(){
450 RealType totalFZ;
451 Vector3d force(0.0);
452 Vector3d com;
453 RealType totalFZ_local = 0;
454 RealType lrPot;
455 std::list<ZconstraintMol>::iterator i;
456 StuntDouble* sd;
457 Molecule::IntegrableObjectIterator ii;
458 Molecule* mol;
459 for ( i = movingZMols_.begin(); i != movingZMols_.end(); ++i) {
460 mol = i->mol;
461 com = mol->getCom();
462 RealType resPos = usingSMD_? i->cantPos : i->param.zTargetPos;
463 RealType diff = com[whichDirection] - resPos;
464 RealType harmonicU = 0.5 * i->param.kz * diff * diff;
465 lrPot = currSnapshot_->getLongRangePotential();
466 lrPot += harmonicU;
467 currSnapshot_->setLongRangePotential(lrPot);
468 RealType harmonicF = -i->param.kz * diff;
469 totalFZ_local += harmonicF;
470
471 //adjust force
472 for(sd = mol->beginIntegrableObject(ii); sd != NULL;
473 sd = mol->nextIntegrableObject(ii)) {
474
475 force[whichDirection] = getHFOfFixedZMols(mol, sd, harmonicF);
476 sd->addFrc(force);
477 }
478 }
479
480 #ifndef IS_MPI
481 totalFZ = totalFZ_local;
482 #else
483 MPI_Allreduce(&totalFZ_local, &totalFZ, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
484 #endif
485
486 //modify the forces of unconstrained molecules
487 std::vector<Molecule*>::iterator j;
488 for ( j = unzconsMols_.begin(); j != unzconsMols_.end(); ++j) {
489
490 mol = *j;
491
492 for(sd = mol->beginIntegrableObject(ii); sd != NULL;
493 sd = mol->nextIntegrableObject(ii)) {
494
495 force[whichDirection] = getHFOfUnconsMols(mol, totalFZ);
496 sd->addFrc(force);
497 }
498 }
499
500 }
501
502 bool ZconstraintForceManager::checkZConsState(){
503 Vector3d com;
504 RealType diff;
505 int changed_local = 0;
506
507 std::list<ZconstraintMol>::iterator i;
508 std::list<ZconstraintMol>::iterator j;
509
510 std::list<ZconstraintMol> newMovingZMols;
511 for ( i = fixedZMols_.begin(); i != fixedZMols_.end();) {
512 com = i->mol->getCom();
513 diff = fabs(com[whichDirection] - i->param.zTargetPos);
514 if (diff > zconsTol_) {
515 if (usingZconsGap_) {
516 i->endFixingTime = infiniteTime;
517 }
518 j = i++;
519 newMovingZMols.push_back(*j);
520 fixedZMols_.erase(j);
521
522 changed_local = 1;
523 }else {
524 ++i;
525 }
526 }
527
528 std::list<ZconstraintMol> newFixedZMols;
529 for ( i = movingZMols_.begin(); i != movingZMols_.end();) {
530 com = i->mol->getCom();
531 diff = fabs(com[whichDirection] - i->param.zTargetPos);
532 if (diff <= zconsTol_) {
533 if (usingZconsGap_) {
534 i->endFixingTime = currSnapshot_->getTime() + zconsFixingTime_;
535 }
536 //this moving zconstraint molecule is about to fixed
537 //moved this molecule to
538 j = i++;
539 newFixedZMols.push_back(*j);
540 movingZMols_.erase(j);
541 changed_local = 1;
542 }else {
543 ++i;
544 }
545 }
546
547 //merge the lists
548 fixedZMols_.insert(fixedZMols_.end(), newFixedZMols.begin(), newFixedZMols.end());
549 movingZMols_.insert(movingZMols_.end(), newMovingZMols.begin(), newMovingZMols.end());
550
551 int changed;
552 #ifndef IS_MPI
553 changed = changed_local;
554 #else
555 MPI_Allreduce(&changed_local, &changed, 1, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
556 #endif
557
558 return (changed > 0);
559 }
560
561 bool ZconstraintForceManager::haveFixedZMols(){
562 int havingFixed;
563 int havingFixed_local = fixedZMols_.empty() ? 0 : 1;
564
565 #ifndef IS_MPI
566 havingFixed = havingFixed_local;
567 #else
568 MPI_Allreduce(&havingFixed_local, &havingFixed, 1, MPI_INT, MPI_SUM,
569 MPI_COMM_WORLD);
570 #endif
571
572 return havingFixed > 0;
573 }
574
575
576 bool ZconstraintForceManager::haveMovingZMols(){
577 int havingMoving_local;
578 int havingMoving;
579
580 havingMoving_local = movingZMols_.empty()? 0 : 1;
581
582 #ifndef IS_MPI
583 havingMoving = havingMoving_local;
584 #else
585 MPI_Allreduce(&havingMoving_local, &havingMoving, 1, MPI_INT, MPI_SUM,
586 MPI_COMM_WORLD);
587 #endif
588
589 return havingMoving > 0;
590 }
591
592 void ZconstraintForceManager::calcTotalMassMovingZMols(){
593
594 RealType totMassMovingZMols_local = 0.0;
595 std::list<ZconstraintMol>::iterator i;
596 for ( i = movingZMols_.begin(); i != movingZMols_.end(); ++i) {
597 totMassMovingZMols_local += i->mol->getMass();
598 }
599
600 #ifdef IS_MPI
601 MPI_Allreduce(&totMassMovingZMols_local, &totMassMovingZMols_, 1, MPI_REALTYPE,
602 MPI_SUM, MPI_COMM_WORLD);
603 #else
604 totMassMovingZMols_ = totMassMovingZMols_local;
605 #endif
606
607 }
608
609 RealType ZconstraintForceManager::getZFOfFixedZMols(Molecule* mol, StuntDouble* sd, RealType totalForce){
610 return totalForce * sd->getMass() / mol->getMass();
611 }
612
613 RealType ZconstraintForceManager::getZFOfMovingMols(Molecule* mol, RealType totalForce){
614 return totalForce * mol->getMass() / (totMassUnconsMols_ + totMassMovingZMols_);
615 }
616
617 RealType ZconstraintForceManager::getHFOfFixedZMols(Molecule* mol, StuntDouble*sd, RealType totalForce){
618 return totalForce * sd->getMass() / mol->getMass();
619 }
620
621 RealType ZconstraintForceManager::getHFOfUnconsMols(Molecule* mol, RealType totalForce){
622 return totalForce * mol->getMass() / totMassUnconsMols_;
623 }
624
625 void ZconstraintForceManager::updateZPos(){
626 RealType curTime = currSnapshot_->getTime();
627 std::list<ZconstraintMol>::iterator i;
628 for ( i = fixedZMols_.begin(); i != fixedZMols_.end(); ++i) {
629 i->param.zTargetPos += zconsGap_;
630 }
631 }
632
633 void ZconstraintForceManager::updateCantPos(){
634 std::list<ZconstraintMol>::iterator i;
635 for ( i = movingZMols_.begin(); i != movingZMols_.end(); ++i) {
636 i->cantPos += i->param.cantVel * dt_;
637 }
638 }
639
640 RealType ZconstraintForceManager::getZTargetPos(int index){
641 RealType zTargetPos;
642 #ifndef IS_MPI
643 Molecule* mol = info_->getMoleculeByGlobalIndex(index);
644 assert(mol);
645 Vector3d com = mol->getCom();
646 zTargetPos = com[whichDirection];
647 #else
648 int whicProc = info_->getMolToProc(index);
649 MPI_Bcast(&zTargetPos, 1, MPI_REALTYPE, whicProc, MPI_COMM_WORLD);
650 #endif
651 return zTargetPos;
652 }
653
654 }

Properties

Name Value
svn:executable *
svn:keywords Author Id Revision Date