ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/constraints/ZconstraintForceManager.cpp
Revision: 1665
Committed: Tue Nov 22 20:38:56 2011 UTC (13 years, 5 months ago) by gezelter
File size: 19642 byte(s)
Log Message:
updated copyright notices

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

Properties

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