ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/constraints/ZconstraintForceManager.cpp
Revision: 2204
Committed: Fri Apr 15 22:04:00 2005 UTC (19 years, 2 months ago) by gezelter
File size: 19468 byte(s)
Log Message:
xemacs has been drafted to perform our indentation services

File Contents

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

Properties

Name Value
svn:executable *