ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/constraints/ZconstraintForceManager.cpp
Revision: 1879
Committed: Sun Jun 16 15:15:42 2013 UTC (11 years, 10 months ago) by gezelter
File size: 19304 byte(s)
Log Message:
MERGE OpenMD development 1783:1878 into trunk

File Contents

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

Properties

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