ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/brains/SimInfo.cpp
Revision: 2419
Committed: Tue Nov 8 13:32:06 2005 UTC (18 years, 8 months ago) by chrisfen
File size: 36730 byte(s)
Log Message:
Added a keyword for output of forces and torques

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     /**
43     * @file SimInfo.cpp
44     * @author tlin
45     * @date 11/02/2004
46     * @version 1.0
47     */
48 gezelter 1490
49 gezelter 1930 #include <algorithm>
50     #include <set>
51 gezelter 1490
52 tim 1492 #include "brains/SimInfo.hpp"
53 gezelter 1930 #include "math/Vector3.hpp"
54     #include "primitives/Molecule.hpp"
55 gezelter 2285 #include "UseTheForce/fCutoffPolicy.h"
56 chrisfen 2305 #include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h"
57 chrisfen 2415 #include "UseTheForce/DarkSide/fElectrostaticScreeningMethod.h"
58 gezelter 1930 #include "UseTheForce/doForces_interface.h"
59 chrisfen 2309 #include "UseTheForce/DarkSide/electrostatic_interface.h"
60 gezelter 1930 #include "UseTheForce/notifyCutoffs_interface.h"
61     #include "utils/MemoryUtils.hpp"
62 tim 1492 #include "utils/simError.h"
63 tim 2000 #include "selection/SelectionManager.hpp"
64 gezelter 1490
65 gezelter 1930 #ifdef IS_MPI
66     #include "UseTheForce/mpiComponentPlan.h"
67     #include "UseTheForce/DarkSide/simParallel_interface.h"
68     #endif
69 gezelter 1490
70 gezelter 1930 namespace oopse {
71 gezelter 1490
72 gezelter 2204 SimInfo::SimInfo(MakeStamps* stamps, std::vector<std::pair<MoleculeStamp*, int> >& molStampPairs,
73     ForceField* ff, Globals* simParams) :
74     stamps_(stamps), forceField_(ff), simParams_(simParams),
75     ndf_(0), ndfRaw_(0), ndfTrans_(0), nZconstraint_(0),
76     nGlobalMols_(0), nGlobalAtoms_(0), nGlobalCutoffGroups_(0),
77     nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0),
78     nAtoms_(0), nBonds_(0), nBends_(0), nTorsions_(0), nRigidBodies_(0),
79     nIntegrableObjects_(0), nCutoffGroups_(0), nConstraints_(0),
80     sman_(NULL), fortranInitialized_(false) {
81 gezelter 1490
82 gezelter 1930
83 gezelter 2204 std::vector<std::pair<MoleculeStamp*, int> >::iterator i;
84     MoleculeStamp* molStamp;
85     int nMolWithSameStamp;
86     int nCutoffAtoms = 0; // number of atoms belong to cutoff groups
87 chrisfen 2344 int nGroups = 0; //total cutoff groups defined in meta-data file
88 gezelter 2204 CutoffGroupStamp* cgStamp;
89     RigidBodyStamp* rbStamp;
90     int nRigidAtoms = 0;
91 gezelter 1930
92 gezelter 2204 for (i = molStampPairs.begin(); i !=molStampPairs.end(); ++i) {
93 gezelter 1930 molStamp = i->first;
94     nMolWithSameStamp = i->second;
95    
96     addMoleculeStamp(molStamp, nMolWithSameStamp);
97 gezelter 1490
98 gezelter 1930 //calculate atoms in molecules
99     nGlobalAtoms_ += molStamp->getNAtoms() *nMolWithSameStamp;
100 gezelter 1490
101    
102 gezelter 1930 //calculate atoms in cutoff groups
103     int nAtomsInGroups = 0;
104     int nCutoffGroupsInStamp = molStamp->getNCutoffGroups();
105    
106     for (int j=0; j < nCutoffGroupsInStamp; j++) {
107 gezelter 2204 cgStamp = molStamp->getCutoffGroup(j);
108     nAtomsInGroups += cgStamp->getNMembers();
109 gezelter 1930 }
110 gezelter 1490
111 gezelter 1930 nGroups += nCutoffGroupsInStamp * nMolWithSameStamp;
112 chrisfen 2344
113 gezelter 1930 nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp;
114 gezelter 1490
115 gezelter 1930 //calculate atoms in rigid bodies
116     int nAtomsInRigidBodies = 0;
117 tim 1958 int nRigidBodiesInStamp = molStamp->getNRigidBodies();
118 gezelter 1930
119     for (int j=0; j < nRigidBodiesInStamp; j++) {
120 gezelter 2204 rbStamp = molStamp->getRigidBody(j);
121     nAtomsInRigidBodies += rbStamp->getNMembers();
122 gezelter 1930 }
123 gezelter 1490
124 gezelter 1930 nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp;
125     nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp;
126    
127 gezelter 2204 }
128 chrisfen 1636
129 chrisfen 2344 //every free atom (atom does not belong to cutoff groups) is a cutoff
130     //group therefore the total number of cutoff groups in the system is
131     //equal to the total number of atoms minus number of atoms belong to
132     //cutoff group defined in meta-data file plus the number of cutoff
133     //groups defined in meta-data file
134 gezelter 2204 nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
135 gezelter 1490
136 chrisfen 2344 //every free atom (atom does not belong to rigid bodies) is an
137     //integrable object therefore the total number of integrable objects
138     //in the system is equal to the total number of atoms minus number of
139     //atoms belong to rigid body defined in meta-data file plus the number
140     //of rigid bodies defined in meta-data file
141     nGlobalIntegrableObjects_ = nGlobalAtoms_ - nRigidAtoms
142     + nGlobalRigidBodies_;
143    
144 gezelter 2204 nGlobalMols_ = molStampIds_.size();
145 gezelter 1490
146 gezelter 1930 #ifdef IS_MPI
147 gezelter 2204 molToProcMap_.resize(nGlobalMols_);
148 gezelter 1930 #endif
149 tim 1976
150 gezelter 2204 }
151 gezelter 1490
152 gezelter 2204 SimInfo::~SimInfo() {
153 tim 2082 std::map<int, Molecule*>::iterator i;
154     for (i = molecules_.begin(); i != molecules_.end(); ++i) {
155 gezelter 2204 delete i->second;
156 tim 2082 }
157     molecules_.clear();
158 tim 2187
159     delete stamps_;
160 gezelter 1930 delete sman_;
161     delete simParams_;
162     delete forceField_;
163 gezelter 2204 }
164 gezelter 1490
165 gezelter 2204 int SimInfo::getNGlobalConstraints() {
166 gezelter 1930 int nGlobalConstraints;
167     #ifdef IS_MPI
168     MPI_Allreduce(&nConstraints_, &nGlobalConstraints, 1, MPI_INT, MPI_SUM,
169     MPI_COMM_WORLD);
170     #else
171     nGlobalConstraints = nConstraints_;
172     #endif
173     return nGlobalConstraints;
174 gezelter 2204 }
175 gezelter 1490
176 gezelter 2204 bool SimInfo::addMolecule(Molecule* mol) {
177 gezelter 1930 MoleculeIterator i;
178 gezelter 1490
179 gezelter 1930 i = molecules_.find(mol->getGlobalIndex());
180     if (i == molecules_.end() ) {
181 gezelter 1490
182 gezelter 2204 molecules_.insert(std::make_pair(mol->getGlobalIndex(), mol));
183 gezelter 1930
184 gezelter 2204 nAtoms_ += mol->getNAtoms();
185     nBonds_ += mol->getNBonds();
186     nBends_ += mol->getNBends();
187     nTorsions_ += mol->getNTorsions();
188     nRigidBodies_ += mol->getNRigidBodies();
189     nIntegrableObjects_ += mol->getNIntegrableObjects();
190     nCutoffGroups_ += mol->getNCutoffGroups();
191     nConstraints_ += mol->getNConstraintPairs();
192 gezelter 1490
193 gezelter 2204 addExcludePairs(mol);
194 gezelter 1930
195 gezelter 2204 return true;
196 gezelter 1930 } else {
197 gezelter 2204 return false;
198 gezelter 1930 }
199 gezelter 2204 }
200 gezelter 1490
201 gezelter 2204 bool SimInfo::removeMolecule(Molecule* mol) {
202 gezelter 1930 MoleculeIterator i;
203     i = molecules_.find(mol->getGlobalIndex());
204 gezelter 1490
205 gezelter 1930 if (i != molecules_.end() ) {
206 gezelter 1490
207 gezelter 2204 assert(mol == i->second);
208 gezelter 1930
209 gezelter 2204 nAtoms_ -= mol->getNAtoms();
210     nBonds_ -= mol->getNBonds();
211     nBends_ -= mol->getNBends();
212     nTorsions_ -= mol->getNTorsions();
213     nRigidBodies_ -= mol->getNRigidBodies();
214     nIntegrableObjects_ -= mol->getNIntegrableObjects();
215     nCutoffGroups_ -= mol->getNCutoffGroups();
216     nConstraints_ -= mol->getNConstraintPairs();
217 gezelter 1490
218 gezelter 2204 removeExcludePairs(mol);
219     molecules_.erase(mol->getGlobalIndex());
220 gezelter 1490
221 gezelter 2204 delete mol;
222 gezelter 1930
223 gezelter 2204 return true;
224 gezelter 1930 } else {
225 gezelter 2204 return false;
226 gezelter 1930 }
227    
228    
229 gezelter 2204 }
230 gezelter 1930
231    
232 gezelter 2204 Molecule* SimInfo::beginMolecule(MoleculeIterator& i) {
233 gezelter 1930 i = molecules_.begin();
234     return i == molecules_.end() ? NULL : i->second;
235 gezelter 2204 }
236 gezelter 1930
237 gezelter 2204 Molecule* SimInfo::nextMolecule(MoleculeIterator& i) {
238 gezelter 1930 ++i;
239     return i == molecules_.end() ? NULL : i->second;
240 gezelter 2204 }
241 gezelter 1490
242    
243 gezelter 2204 void SimInfo::calcNdf() {
244 gezelter 1930 int ndf_local;
245     MoleculeIterator i;
246     std::vector<StuntDouble*>::iterator j;
247     Molecule* mol;
248     StuntDouble* integrableObject;
249 gezelter 1490
250 gezelter 1930 ndf_local = 0;
251    
252     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
253 gezelter 2204 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
254     integrableObject = mol->nextIntegrableObject(j)) {
255 gezelter 1490
256 gezelter 2204 ndf_local += 3;
257 gezelter 1490
258 gezelter 2204 if (integrableObject->isDirectional()) {
259     if (integrableObject->isLinear()) {
260     ndf_local += 2;
261     } else {
262     ndf_local += 3;
263     }
264     }
265 gezelter 1930
266 gezelter 2204 }//end for (integrableObject)
267 gezelter 1930 }// end for (mol)
268    
269     // n_constraints is local, so subtract them on each processor
270     ndf_local -= nConstraints_;
271    
272     #ifdef IS_MPI
273     MPI_Allreduce(&ndf_local,&ndf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
274     #else
275     ndf_ = ndf_local;
276     #endif
277    
278     // nZconstraints_ is global, as are the 3 COM translations for the
279     // entire system:
280     ndf_ = ndf_ - 3 - nZconstraint_;
281    
282 gezelter 2204 }
283 gezelter 1490
284 gezelter 2204 void SimInfo::calcNdfRaw() {
285 gezelter 1930 int ndfRaw_local;
286 gezelter 1490
287 gezelter 1930 MoleculeIterator i;
288     std::vector<StuntDouble*>::iterator j;
289     Molecule* mol;
290     StuntDouble* integrableObject;
291    
292     // Raw degrees of freedom that we have to set
293     ndfRaw_local = 0;
294    
295     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
296 gezelter 2204 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
297     integrableObject = mol->nextIntegrableObject(j)) {
298 gezelter 1930
299 gezelter 2204 ndfRaw_local += 3;
300 gezelter 1930
301 gezelter 2204 if (integrableObject->isDirectional()) {
302     if (integrableObject->isLinear()) {
303     ndfRaw_local += 2;
304     } else {
305     ndfRaw_local += 3;
306     }
307     }
308 gezelter 1930
309 gezelter 2204 }
310 gezelter 1930 }
311    
312     #ifdef IS_MPI
313     MPI_Allreduce(&ndfRaw_local,&ndfRaw_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
314     #else
315     ndfRaw_ = ndfRaw_local;
316     #endif
317 gezelter 2204 }
318 gezelter 1490
319 gezelter 2204 void SimInfo::calcNdfTrans() {
320 gezelter 1930 int ndfTrans_local;
321 gezelter 1490
322 gezelter 1930 ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_;
323 gezelter 1490
324    
325 gezelter 1930 #ifdef IS_MPI
326     MPI_Allreduce(&ndfTrans_local,&ndfTrans_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
327     #else
328     ndfTrans_ = ndfTrans_local;
329     #endif
330 gezelter 1490
331 gezelter 1930 ndfTrans_ = ndfTrans_ - 3 - nZconstraint_;
332    
333 gezelter 2204 }
334 gezelter 1490
335 gezelter 2204 void SimInfo::addExcludePairs(Molecule* mol) {
336 gezelter 1930 std::vector<Bond*>::iterator bondIter;
337     std::vector<Bend*>::iterator bendIter;
338     std::vector<Torsion*>::iterator torsionIter;
339     Bond* bond;
340     Bend* bend;
341     Torsion* torsion;
342     int a;
343     int b;
344     int c;
345     int d;
346    
347     for (bond= mol->beginBond(bondIter); bond != NULL; bond = mol->nextBond(bondIter)) {
348 gezelter 2204 a = bond->getAtomA()->getGlobalIndex();
349     b = bond->getAtomB()->getGlobalIndex();
350     exclude_.addPair(a, b);
351 gezelter 1930 }
352 gezelter 1490
353 gezelter 1930 for (bend= mol->beginBend(bendIter); bend != NULL; bend = mol->nextBend(bendIter)) {
354 gezelter 2204 a = bend->getAtomA()->getGlobalIndex();
355     b = bend->getAtomB()->getGlobalIndex();
356     c = bend->getAtomC()->getGlobalIndex();
357 gezelter 1490
358 gezelter 2204 exclude_.addPair(a, b);
359     exclude_.addPair(a, c);
360     exclude_.addPair(b, c);
361 gezelter 1930 }
362 gezelter 1490
363 gezelter 1930 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL; torsion = mol->nextTorsion(torsionIter)) {
364 gezelter 2204 a = torsion->getAtomA()->getGlobalIndex();
365     b = torsion->getAtomB()->getGlobalIndex();
366     c = torsion->getAtomC()->getGlobalIndex();
367     d = torsion->getAtomD()->getGlobalIndex();
368 gezelter 1490
369 gezelter 2204 exclude_.addPair(a, b);
370     exclude_.addPair(a, c);
371     exclude_.addPair(a, d);
372     exclude_.addPair(b, c);
373     exclude_.addPair(b, d);
374     exclude_.addPair(c, d);
375 gezelter 1490 }
376    
377 tim 2114 Molecule::RigidBodyIterator rbIter;
378     RigidBody* rb;
379     for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
380 gezelter 2204 std::vector<Atom*> atoms = rb->getAtoms();
381     for (int i = 0; i < atoms.size() -1 ; ++i) {
382     for (int j = i + 1; j < atoms.size(); ++j) {
383     a = atoms[i]->getGlobalIndex();
384     b = atoms[j]->getGlobalIndex();
385     exclude_.addPair(a, b);
386     }
387     }
388 tim 2114 }
389    
390 gezelter 2204 }
391 gezelter 1930
392 gezelter 2204 void SimInfo::removeExcludePairs(Molecule* mol) {
393 gezelter 1930 std::vector<Bond*>::iterator bondIter;
394     std::vector<Bend*>::iterator bendIter;
395     std::vector<Torsion*>::iterator torsionIter;
396     Bond* bond;
397     Bend* bend;
398     Torsion* torsion;
399     int a;
400     int b;
401     int c;
402     int d;
403    
404     for (bond= mol->beginBond(bondIter); bond != NULL; bond = mol->nextBond(bondIter)) {
405 gezelter 2204 a = bond->getAtomA()->getGlobalIndex();
406     b = bond->getAtomB()->getGlobalIndex();
407     exclude_.removePair(a, b);
408 gezelter 1490 }
409 gezelter 1930
410     for (bend= mol->beginBend(bendIter); bend != NULL; bend = mol->nextBend(bendIter)) {
411 gezelter 2204 a = bend->getAtomA()->getGlobalIndex();
412     b = bend->getAtomB()->getGlobalIndex();
413     c = bend->getAtomC()->getGlobalIndex();
414 gezelter 1930
415 gezelter 2204 exclude_.removePair(a, b);
416     exclude_.removePair(a, c);
417     exclude_.removePair(b, c);
418 gezelter 1490 }
419 gezelter 1930
420     for (torsion= mol->beginTorsion(torsionIter); torsion != NULL; torsion = mol->nextTorsion(torsionIter)) {
421 gezelter 2204 a = torsion->getAtomA()->getGlobalIndex();
422     b = torsion->getAtomB()->getGlobalIndex();
423     c = torsion->getAtomC()->getGlobalIndex();
424     d = torsion->getAtomD()->getGlobalIndex();
425 gezelter 1930
426 gezelter 2204 exclude_.removePair(a, b);
427     exclude_.removePair(a, c);
428     exclude_.removePair(a, d);
429     exclude_.removePair(b, c);
430     exclude_.removePair(b, d);
431     exclude_.removePair(c, d);
432 gezelter 1930 }
433    
434 tim 2114 Molecule::RigidBodyIterator rbIter;
435     RigidBody* rb;
436     for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
437 gezelter 2204 std::vector<Atom*> atoms = rb->getAtoms();
438     for (int i = 0; i < atoms.size() -1 ; ++i) {
439     for (int j = i + 1; j < atoms.size(); ++j) {
440     a = atoms[i]->getGlobalIndex();
441     b = atoms[j]->getGlobalIndex();
442     exclude_.removePair(a, b);
443     }
444     }
445 tim 2114 }
446    
447 gezelter 2204 }
448 gezelter 1490
449    
450 gezelter 2204 void SimInfo::addMoleculeStamp(MoleculeStamp* molStamp, int nmol) {
451 gezelter 1930 int curStampId;
452 gezelter 1490
453 gezelter 1930 //index from 0
454     curStampId = moleculeStamps_.size();
455 gezelter 1490
456 gezelter 1930 moleculeStamps_.push_back(molStamp);
457     molStampIds_.insert(molStampIds_.end(), nmol, curStampId);
458 gezelter 2204 }
459 gezelter 1490
460 gezelter 2204 void SimInfo::update() {
461 gezelter 1490
462 gezelter 1930 setupSimType();
463 gezelter 1490
464 gezelter 1930 #ifdef IS_MPI
465     setupFortranParallel();
466     #endif
467 gezelter 1490
468 gezelter 1930 setupFortranSim();
469 gezelter 1490
470 gezelter 1930 //setup fortran force field
471     /** @deprecate */
472     int isError = 0;
473 chrisfen 2297
474 chrisfen 2302 setupElectrostaticSummationMethod( isError );
475 chrisfen 2297
476 gezelter 1930 if(isError){
477 gezelter 2204 sprintf( painCave.errMsg,
478     "ForceField error: There was an error initializing the forceField in fortran.\n" );
479     painCave.isFatal = 1;
480     simError();
481 gezelter 1930 }
482 gezelter 1490
483 gezelter 1930
484     setupCutoff();
485 gezelter 1490
486 gezelter 1930 calcNdf();
487     calcNdfRaw();
488     calcNdfTrans();
489    
490     fortranInitialized_ = true;
491 gezelter 2204 }
492 gezelter 1490
493 gezelter 2204 std::set<AtomType*> SimInfo::getUniqueAtomTypes() {
494 gezelter 1930 SimInfo::MoleculeIterator mi;
495     Molecule* mol;
496     Molecule::AtomIterator ai;
497     Atom* atom;
498     std::set<AtomType*> atomTypes;
499 gezelter 1490
500 gezelter 1930 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
501 gezelter 1490
502 gezelter 2204 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
503     atomTypes.insert(atom->getAtomType());
504     }
505 gezelter 1930
506     }
507 gezelter 1490
508 gezelter 1930 return atomTypes;
509 gezelter 2204 }
510 gezelter 1490
511 gezelter 2204 void SimInfo::setupSimType() {
512 gezelter 1930 std::set<AtomType*>::iterator i;
513     std::set<AtomType*> atomTypes;
514     atomTypes = getUniqueAtomTypes();
515 gezelter 1490
516 gezelter 1930 int useLennardJones = 0;
517     int useElectrostatic = 0;
518     int useEAM = 0;
519     int useCharge = 0;
520     int useDirectional = 0;
521     int useDipole = 0;
522     int useGayBerne = 0;
523     int useSticky = 0;
524 chrisfen 2220 int useStickyPower = 0;
525 gezelter 1930 int useShape = 0;
526     int useFLARB = 0; //it is not in AtomType yet
527     int useDirectionalAtom = 0;
528     int useElectrostatics = 0;
529     //usePBC and useRF are from simParams
530 tim 2364 int usePBC = simParams_->getUsePeriodicBoundaryConditions();
531 chrisfen 2310 int useRF;
532 chrisfen 2419 int useSF;
533 tim 2364 std::string myMethod;
534 gezelter 1490
535 chrisfen 2310 // set the useRF logical
536 tim 2364 useRF = 0;
537 chrisfen 2419 useSF = 0;
538 chrisfen 2390
539    
540 tim 2364 if (simParams_->haveElectrostaticSummationMethod()) {
541 chrisfen 2390 std::string myMethod = simParams_->getElectrostaticSummationMethod();
542     toUpper(myMethod);
543     if (myMethod == "REACTION_FIELD") {
544     useRF=1;
545 chrisfen 2404 } else {
546 chrisfen 2419 if (myMethod == "SHIFTED_FORCE") {
547     useSF = 1;
548 chrisfen 2404 }
549 chrisfen 2390 }
550 tim 2364 }
551 chrisfen 2310
552 gezelter 1930 //loop over all of the atom types
553     for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
554 gezelter 2204 useLennardJones |= (*i)->isLennardJones();
555     useElectrostatic |= (*i)->isElectrostatic();
556     useEAM |= (*i)->isEAM();
557     useCharge |= (*i)->isCharge();
558     useDirectional |= (*i)->isDirectional();
559     useDipole |= (*i)->isDipole();
560     useGayBerne |= (*i)->isGayBerne();
561     useSticky |= (*i)->isSticky();
562 chrisfen 2220 useStickyPower |= (*i)->isStickyPower();
563 gezelter 2204 useShape |= (*i)->isShape();
564 gezelter 1930 }
565 gezelter 1490
566 chrisfen 2220 if (useSticky || useStickyPower || useDipole || useGayBerne || useShape) {
567 gezelter 2204 useDirectionalAtom = 1;
568 gezelter 1930 }
569 gezelter 1490
570 gezelter 1930 if (useCharge || useDipole) {
571 gezelter 2204 useElectrostatics = 1;
572 gezelter 1930 }
573 gezelter 1490
574 gezelter 1930 #ifdef IS_MPI
575     int temp;
576 gezelter 1490
577 gezelter 1930 temp = usePBC;
578     MPI_Allreduce(&temp, &usePBC, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
579 gezelter 1490
580 gezelter 1930 temp = useDirectionalAtom;
581     MPI_Allreduce(&temp, &useDirectionalAtom, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
582 gezelter 1490
583 gezelter 1930 temp = useLennardJones;
584     MPI_Allreduce(&temp, &useLennardJones, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
585 gezelter 1490
586 gezelter 1930 temp = useElectrostatics;
587     MPI_Allreduce(&temp, &useElectrostatics, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
588 gezelter 1490
589 gezelter 1930 temp = useCharge;
590     MPI_Allreduce(&temp, &useCharge, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
591 gezelter 1490
592 gezelter 1930 temp = useDipole;
593     MPI_Allreduce(&temp, &useDipole, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
594 gezelter 1490
595 gezelter 1930 temp = useSticky;
596     MPI_Allreduce(&temp, &useSticky, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
597 gezelter 1490
598 chrisfen 2220 temp = useStickyPower;
599     MPI_Allreduce(&temp, &useStickyPower, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
600    
601 gezelter 1930 temp = useGayBerne;
602     MPI_Allreduce(&temp, &useGayBerne, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
603 gezelter 1490
604 gezelter 1930 temp = useEAM;
605     MPI_Allreduce(&temp, &useEAM, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
606 gezelter 1490
607 gezelter 1930 temp = useShape;
608     MPI_Allreduce(&temp, &useShape, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
609    
610     temp = useFLARB;
611     MPI_Allreduce(&temp, &useFLARB, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
612    
613 chrisfen 2310 temp = useRF;
614     MPI_Allreduce(&temp, &useRF, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
615    
616 chrisfen 2419 temp = useSF;
617     MPI_Allreduce(&temp, &useSF, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
618 chrisfen 2404
619 gezelter 1490 #endif
620    
621 gezelter 1930 fInfo_.SIM_uses_PBC = usePBC;
622     fInfo_.SIM_uses_DirectionalAtoms = useDirectionalAtom;
623     fInfo_.SIM_uses_LennardJones = useLennardJones;
624     fInfo_.SIM_uses_Electrostatics = useElectrostatics;
625     fInfo_.SIM_uses_Charges = useCharge;
626     fInfo_.SIM_uses_Dipoles = useDipole;
627     fInfo_.SIM_uses_Sticky = useSticky;
628 chrisfen 2220 fInfo_.SIM_uses_StickyPower = useStickyPower;
629 gezelter 1930 fInfo_.SIM_uses_GayBerne = useGayBerne;
630     fInfo_.SIM_uses_EAM = useEAM;
631     fInfo_.SIM_uses_Shapes = useShape;
632     fInfo_.SIM_uses_FLARB = useFLARB;
633 chrisfen 2310 fInfo_.SIM_uses_RF = useRF;
634 chrisfen 2419 fInfo_.SIM_uses_SF = useSF;
635 gezelter 1490
636 chrisfen 2390 if( myMethod == "REACTION_FIELD") {
637    
638 gezelter 2204 if (simParams_->haveDielectric()) {
639     fInfo_.dielect = simParams_->getDielectric();
640     } else {
641     sprintf(painCave.errMsg,
642     "SimSetup Error: No Dielectric constant was set.\n"
643     "\tYou are trying to use Reaction Field without"
644     "\tsetting a dielectric constant!\n");
645     painCave.isFatal = 1;
646     simError();
647 chrisfen 2390 }
648 gezelter 1930 }
649 chrisfen 2404
650 gezelter 2204 }
651 gezelter 1490
652 gezelter 2204 void SimInfo::setupFortranSim() {
653 gezelter 1930 int isError;
654     int nExclude;
655     std::vector<int> fortranGlobalGroupMembership;
656    
657     nExclude = exclude_.getSize();
658     isError = 0;
659 gezelter 1490
660 gezelter 1930 //globalGroupMembership_ is filled by SimCreator
661     for (int i = 0; i < nGlobalAtoms_; i++) {
662 gezelter 2204 fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1);
663 gezelter 1930 }
664 gezelter 1490
665 gezelter 1930 //calculate mass ratio of cutoff group
666     std::vector<double> mfact;
667     SimInfo::MoleculeIterator mi;
668     Molecule* mol;
669     Molecule::CutoffGroupIterator ci;
670     CutoffGroup* cg;
671     Molecule::AtomIterator ai;
672     Atom* atom;
673     double totalMass;
674    
675     //to avoid memory reallocation, reserve enough space for mfact
676     mfact.reserve(getNCutoffGroups());
677 gezelter 1490
678 gezelter 1930 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
679 gezelter 2204 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
680 gezelter 1490
681 gezelter 2204 totalMass = cg->getMass();
682     for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
683 chrisfen 2344 // Check for massless groups - set mfact to 1 if true
684     if (totalMass != 0)
685     mfact.push_back(atom->getMass()/totalMass);
686     else
687     mfact.push_back( 1.0 );
688 gezelter 2204 }
689 gezelter 1490
690 gezelter 2204 }
691 gezelter 1930 }
692 gezelter 1490
693 gezelter 1930 //fill ident array of local atoms (it is actually ident of AtomType, it is so confusing !!!)
694     std::vector<int> identArray;
695 gezelter 1490
696 gezelter 1930 //to avoid memory reallocation, reserve enough space identArray
697     identArray.reserve(getNAtoms());
698    
699     for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
700 gezelter 2204 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
701     identArray.push_back(atom->getIdent());
702     }
703 gezelter 1930 }
704 gezelter 1490
705 gezelter 1930 //fill molMembershipArray
706     //molMembershipArray is filled by SimCreator
707     std::vector<int> molMembershipArray(nGlobalAtoms_);
708     for (int i = 0; i < nGlobalAtoms_; i++) {
709 gezelter 2204 molMembershipArray[i] = globalMolMembership_[i] + 1;
710 gezelter 1930 }
711    
712     //setup fortran simulation
713     int nGlobalExcludes = 0;
714     int* globalExcludes = NULL;
715     int* excludeList = exclude_.getExcludeList();
716     setFortranSim( &fInfo_, &nGlobalAtoms_, &nAtoms_, &identArray[0], &nExclude, excludeList ,
717 gezelter 2204 &nGlobalExcludes, globalExcludes, &molMembershipArray[0],
718     &mfact[0], &nCutoffGroups_, &fortranGlobalGroupMembership[0], &isError);
719 gezelter 1490
720 gezelter 1930 if( isError ){
721 gezelter 1490
722 gezelter 2204 sprintf( painCave.errMsg,
723     "There was an error setting the simulation information in fortran.\n" );
724     painCave.isFatal = 1;
725     painCave.severity = OOPSE_ERROR;
726     simError();
727 gezelter 1930 }
728    
729     #ifdef IS_MPI
730     sprintf( checkPointMsg,
731 gezelter 2204 "succesfully sent the simulation information to fortran.\n");
732 gezelter 1930 MPIcheckPoint();
733     #endif // is_mpi
734 gezelter 2204 }
735 gezelter 1490
736    
737 gezelter 1930 #ifdef IS_MPI
738 gezelter 2204 void SimInfo::setupFortranParallel() {
739 gezelter 1930
740     //SimInfo is responsible for creating localToGlobalAtomIndex and localToGlobalGroupIndex
741     std::vector<int> localToGlobalAtomIndex(getNAtoms(), 0);
742     std::vector<int> localToGlobalCutoffGroupIndex;
743     SimInfo::MoleculeIterator mi;
744     Molecule::AtomIterator ai;
745     Molecule::CutoffGroupIterator ci;
746     Molecule* mol;
747     Atom* atom;
748     CutoffGroup* cg;
749     mpiSimData parallelData;
750     int isError;
751 gezelter 1490
752 gezelter 1930 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
753 gezelter 1490
754 gezelter 2204 //local index(index in DataStorge) of atom is important
755     for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
756     localToGlobalAtomIndex[atom->getLocalIndex()] = atom->getGlobalIndex() + 1;
757     }
758 gezelter 1490
759 gezelter 2204 //local index of cutoff group is trivial, it only depends on the order of travesing
760     for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
761     localToGlobalCutoffGroupIndex.push_back(cg->getGlobalIndex() + 1);
762     }
763 gezelter 1930
764     }
765 gezelter 1490
766 gezelter 1930 //fill up mpiSimData struct
767     parallelData.nMolGlobal = getNGlobalMolecules();
768     parallelData.nMolLocal = getNMolecules();
769     parallelData.nAtomsGlobal = getNGlobalAtoms();
770     parallelData.nAtomsLocal = getNAtoms();
771     parallelData.nGroupsGlobal = getNGlobalCutoffGroups();
772     parallelData.nGroupsLocal = getNCutoffGroups();
773     parallelData.myNode = worldRank;
774     MPI_Comm_size(MPI_COMM_WORLD, &(parallelData.nProcessors));
775 gezelter 1490
776 gezelter 1930 //pass mpiSimData struct and index arrays to fortran
777     setFsimParallel(&parallelData, &(parallelData.nAtomsLocal),
778     &localToGlobalAtomIndex[0], &(parallelData.nGroupsLocal),
779     &localToGlobalCutoffGroupIndex[0], &isError);
780 gezelter 1490
781 gezelter 1930 if (isError) {
782 gezelter 2204 sprintf(painCave.errMsg,
783     "mpiRefresh errror: fortran didn't like something we gave it.\n");
784     painCave.isFatal = 1;
785     simError();
786 gezelter 1930 }
787 gezelter 1490
788 gezelter 1930 sprintf(checkPointMsg, " mpiRefresh successful.\n");
789     MPIcheckPoint();
790 gezelter 1490
791    
792 gezelter 2204 }
793 chrisfen 1636
794 gezelter 1930 #endif
795 chrisfen 1636
796 gezelter 2204 double SimInfo::calcMaxCutoffRadius() {
797 chrisfen 1636
798    
799 gezelter 1930 std::set<AtomType*> atomTypes;
800     std::set<AtomType*>::iterator i;
801     std::vector<double> cutoffRadius;
802 gezelter 1490
803 gezelter 1930 //get the unique atom types
804     atomTypes = getUniqueAtomTypes();
805    
806     //query the max cutoff radius among these atom types
807     for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
808 gezelter 2204 cutoffRadius.push_back(forceField_->getRcutFromAtomType(*i));
809 gezelter 1930 }
810    
811     double maxCutoffRadius = *(std::max_element(cutoffRadius.begin(), cutoffRadius.end()));
812 gezelter 1490 #ifdef IS_MPI
813 gezelter 1930 //pick the max cutoff radius among the processors
814 gezelter 1490 #endif
815    
816 gezelter 1930 return maxCutoffRadius;
817 gezelter 2204 }
818 gezelter 1930
819 gezelter 2204 void SimInfo::getCutoff(double& rcut, double& rsw) {
820 gezelter 1490
821 gezelter 1930 if (fInfo_.SIM_uses_Charges | fInfo_.SIM_uses_Dipoles | fInfo_.SIM_uses_RF) {
822    
823 tim 2364 if (!simParams_->haveCutoffRadius()){
824 gezelter 2204 sprintf(painCave.errMsg,
825 gezelter 1930 "SimCreator Warning: No value was set for the cutoffRadius.\n"
826     "\tOOPSE will use a default value of 15.0 angstroms"
827     "\tfor the cutoffRadius.\n");
828 gezelter 2204 painCave.isFatal = 0;
829     simError();
830     rcut = 15.0;
831     } else{
832 tim 2364 rcut = simParams_->getCutoffRadius();
833 gezelter 2204 }
834 gezelter 1930
835 tim 2364 if (!simParams_->haveSwitchingRadius()){
836 gezelter 2204 sprintf(painCave.errMsg,
837 gezelter 1930 "SimCreator Warning: No value was set for switchingRadius.\n"
838     "\tOOPSE will use a default value of\n"
839 chrisfen 2400 "\t0.85 * cutoffRadius for the switchingRadius\n");
840 gezelter 2204 painCave.isFatal = 0;
841     simError();
842 chrisfen 2400 rsw = 0.85 * rcut;
843 gezelter 2204 } else{
844 tim 2364 rsw = simParams_->getSwitchingRadius();
845 gezelter 2204 }
846 gezelter 1930
847     } else {
848 gezelter 2204 // if charge, dipole or reaction field is not used and the cutofff radius is not specified in
849     //meta-data file, the maximum cutoff radius calculated from forcefiled will be used
850 gezelter 1930
851 tim 2364 if (simParams_->haveCutoffRadius()) {
852     rcut = simParams_->getCutoffRadius();
853 gezelter 2204 } else {
854     //set cutoff radius to the maximum cutoff radius based on atom types in the whole system
855     rcut = calcMaxCutoffRadius();
856     }
857 gezelter 1930
858 tim 2364 if (simParams_->haveSwitchingRadius()) {
859     rsw = simParams_->getSwitchingRadius();
860 gezelter 2204 } else {
861     rsw = rcut;
862     }
863 gezelter 1930
864     }
865 gezelter 2204 }
866 tim 2010
867 gezelter 2285 void SimInfo::setupCutoff() {
868 tim 2010 getCutoff(rcut_, rsw_);
869 gezelter 1930 double rnblist = rcut_ + 1; // skin of neighbor list
870    
871     //Pass these cutoff radius etc. to fortran. This function should be called once and only once
872 gezelter 2285
873     int cp = TRADITIONAL_CUTOFF_POLICY;
874     if (simParams_->haveCutoffPolicy()) {
875     std::string myPolicy = simParams_->getCutoffPolicy();
876 tim 2364 toUpper(myPolicy);
877 gezelter 2285 if (myPolicy == "MIX") {
878     cp = MIX_CUTOFF_POLICY;
879     } else {
880     if (myPolicy == "MAX") {
881     cp = MAX_CUTOFF_POLICY;
882     } else {
883     if (myPolicy == "TRADITIONAL") {
884     cp = TRADITIONAL_CUTOFF_POLICY;
885     } else {
886     // throw error
887     sprintf( painCave.errMsg,
888     "SimInfo error: Unknown cutoffPolicy. (Input file specified %s .)\n\tcutoffPolicy must be one of: \"Mix\", \"Max\", or \"Traditional\".", myPolicy.c_str() );
889     painCave.isFatal = 1;
890     simError();
891     }
892     }
893     }
894     }
895 chuckv 2328
896    
897     if (simParams_->haveSkinThickness()) {
898     double skinThickness = simParams_->getSkinThickness();
899     }
900    
901 gezelter 2285 notifyFortranCutoffs(&rcut_, &rsw_, &rnblist, &cp);
902 chrisfen 2309 // also send cutoff notification to electrostatics
903 chrisfen 2381 setElectrostaticCutoffRadius(&rcut_, &rsw_);
904 gezelter 2204 }
905 gezelter 1490
906 chrisfen 2302 void SimInfo::setupElectrostaticSummationMethod( int isError ) {
907 chrisfen 2297
908     int errorOut;
909 chrisfen 2302 int esm = NONE;
910 chrisfen 2408 int sm = UNDAMPED;
911 chrisfen 2297 double alphaVal;
912 chrisfen 2309 double dielectric;
913 chrisfen 2297
914     errorOut = isError;
915 chrisfen 2309 alphaVal = simParams_->getDampingAlpha();
916     dielectric = simParams_->getDielectric();
917 chrisfen 2297
918 chrisfen 2302 if (simParams_->haveElectrostaticSummationMethod()) {
919 chrisfen 2303 std::string myMethod = simParams_->getElectrostaticSummationMethod();
920 tim 2364 toUpper(myMethod);
921 chrisfen 2302 if (myMethod == "NONE") {
922     esm = NONE;
923 chrisfen 2297 } else {
924 chrisfen 2408 if (myMethod == "SWITCHING_FUNCTION") {
925     esm = SWITCHING_FUNCTION;
926 chrisfen 2297 } else {
927 chrisfen 2408 if (myMethod == "SHIFTED_POTENTIAL") {
928     esm = SHIFTED_POTENTIAL;
929     } else {
930     if (myMethod == "SHIFTED_FORCE") {
931     esm = SHIFTED_FORCE;
932 chrisfen 2297 } else {
933 chrisfen 2408 if (myMethod == "REACTION_FIELD") {
934     esm = REACTION_FIELD;
935     } else {
936     // throw error
937     sprintf( painCave.errMsg,
938     "SimInfo error: Unknown electrostaticSummationMethod. (Input file specified %s .)\n\telectrostaticSummationMethod must be one of: \"none\", \"shifted_potential\", \"shifted_force\", or \"reaction_field\".", myMethod.c_str() );
939     painCave.isFatal = 1;
940     simError();
941     }
942     }
943     }
944 chrisfen 2297 }
945     }
946     }
947 chrisfen 2408
948 chrisfen 2415 if (simParams_->haveElectrostaticScreeningMethod()) {
949     std::string myScreen = simParams_->getElectrostaticScreeningMethod();
950 chrisfen 2408 toUpper(myScreen);
951     if (myScreen == "UNDAMPED") {
952     sm = UNDAMPED;
953     } else {
954     if (myScreen == "DAMPED") {
955     sm = DAMPED;
956     if (!simParams_->haveDampingAlpha()) {
957     //throw error
958     sprintf( painCave.errMsg,
959     "SimInfo warning: dampingAlpha was not specified in the input file. A default value of %f (1/ang) will be used.", alphaVal);
960     painCave.isFatal = 0;
961     simError();
962     }
963 chrisfen 2415 } else {
964     // throw error
965     sprintf( painCave.errMsg,
966     "SimInfo error: Unknown electrostaticScreeningMethod. (Input file specified %s .)\n\telectrostaticScreeningMethod must be one of: \"undamped\" or \"damped\".", myScreen.c_str() );
967     painCave.isFatal = 1;
968     simError();
969 chrisfen 2408 }
970     }
971     }
972 chrisfen 2415
973 chrisfen 2309 // let's pass some summation method variables to fortran
974     setElectrostaticSummationMethod( &esm );
975 chrisfen 2408 setScreeningMethod( &sm );
976     setDampingAlpha( &alphaVal );
977 chrisfen 2309 setReactionFieldDielectric( &dielectric );
978     initFortranFF( &esm, &errorOut );
979 chrisfen 2297 }
980    
981 gezelter 2204 void SimInfo::addProperty(GenericData* genData) {
982 gezelter 1930 properties_.addProperty(genData);
983 gezelter 2204 }
984 gezelter 1490
985 gezelter 2204 void SimInfo::removeProperty(const std::string& propName) {
986 gezelter 1930 properties_.removeProperty(propName);
987 gezelter 2204 }
988 gezelter 1490
989 gezelter 2204 void SimInfo::clearProperties() {
990 gezelter 1930 properties_.clearProperties();
991 gezelter 2204 }
992 gezelter 1490
993 gezelter 2204 std::vector<std::string> SimInfo::getPropertyNames() {
994 gezelter 1930 return properties_.getPropertyNames();
995 gezelter 2204 }
996 gezelter 1930
997 gezelter 2204 std::vector<GenericData*> SimInfo::getProperties() {
998 gezelter 1930 return properties_.getProperties();
999 gezelter 2204 }
1000 gezelter 1490
1001 gezelter 2204 GenericData* SimInfo::getPropertyByName(const std::string& propName) {
1002 gezelter 1930 return properties_.getPropertyByName(propName);
1003 gezelter 2204 }
1004 gezelter 1490
1005 gezelter 2204 void SimInfo::setSnapshotManager(SnapshotManager* sman) {
1006 tim 2116 if (sman_ == sman) {
1007 gezelter 2204 return;
1008 tim 2116 }
1009     delete sman_;
1010 gezelter 1930 sman_ = sman;
1011 gezelter 1490
1012 gezelter 1930 Molecule* mol;
1013     RigidBody* rb;
1014     Atom* atom;
1015     SimInfo::MoleculeIterator mi;
1016     Molecule::RigidBodyIterator rbIter;
1017     Molecule::AtomIterator atomIter;;
1018    
1019     for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
1020    
1021 gezelter 2204 for (atom = mol->beginAtom(atomIter); atom != NULL; atom = mol->nextAtom(atomIter)) {
1022     atom->setSnapshotManager(sman_);
1023     }
1024 gezelter 1930
1025 gezelter 2204 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
1026     rb->setSnapshotManager(sman_);
1027     }
1028 gezelter 1930 }
1029 gezelter 1490
1030 gezelter 2204 }
1031 gezelter 1490
1032 gezelter 2204 Vector3d SimInfo::getComVel(){
1033 gezelter 1930 SimInfo::MoleculeIterator i;
1034     Molecule* mol;
1035 gezelter 1490
1036 gezelter 1930 Vector3d comVel(0.0);
1037     double totalMass = 0.0;
1038 gezelter 1490
1039 gezelter 1930
1040     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1041 gezelter 2204 double mass = mol->getMass();
1042     totalMass += mass;
1043     comVel += mass * mol->getComVel();
1044 gezelter 1930 }
1045 gezelter 1490
1046 gezelter 1930 #ifdef IS_MPI
1047     double tmpMass = totalMass;
1048     Vector3d tmpComVel(comVel);
1049     MPI_Allreduce(&tmpMass,&totalMass,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1050     MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1051     #endif
1052    
1053     comVel /= totalMass;
1054    
1055     return comVel;
1056 gezelter 2204 }
1057 gezelter 1490
1058 gezelter 2204 Vector3d SimInfo::getCom(){
1059 gezelter 1930 SimInfo::MoleculeIterator i;
1060     Molecule* mol;
1061 gezelter 1490
1062 gezelter 1930 Vector3d com(0.0);
1063     double totalMass = 0.0;
1064    
1065     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1066 gezelter 2204 double mass = mol->getMass();
1067     totalMass += mass;
1068     com += mass * mol->getCom();
1069 gezelter 1930 }
1070 gezelter 1490
1071     #ifdef IS_MPI
1072 gezelter 1930 double tmpMass = totalMass;
1073     Vector3d tmpCom(com);
1074     MPI_Allreduce(&tmpMass,&totalMass,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1075     MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1076 gezelter 1490 #endif
1077    
1078 gezelter 1930 com /= totalMass;
1079 gezelter 1490
1080 gezelter 1930 return com;
1081 gezelter 1490
1082 gezelter 2204 }
1083 gezelter 1930
1084 gezelter 2204 std::ostream& operator <<(std::ostream& o, SimInfo& info) {
1085 gezelter 1930
1086     return o;
1087 gezelter 2204 }
1088 chuckv 2252
1089    
1090     /*
1091     Returns center of mass and center of mass velocity in one function call.
1092     */
1093    
1094     void SimInfo::getComAll(Vector3d &com, Vector3d &comVel){
1095     SimInfo::MoleculeIterator i;
1096     Molecule* mol;
1097    
1098    
1099     double totalMass = 0.0;
1100    
1101 gezelter 1930
1102 chuckv 2252 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1103     double mass = mol->getMass();
1104     totalMass += mass;
1105     com += mass * mol->getCom();
1106     comVel += mass * mol->getComVel();
1107     }
1108    
1109     #ifdef IS_MPI
1110     double tmpMass = totalMass;
1111     Vector3d tmpCom(com);
1112     Vector3d tmpComVel(comVel);
1113     MPI_Allreduce(&tmpMass,&totalMass,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1114     MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1115     MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1116     #endif
1117    
1118     com /= totalMass;
1119     comVel /= totalMass;
1120     }
1121    
1122     /*
1123     Return intertia tensor for entire system and angular momentum Vector.
1124 chuckv 2256
1125    
1126     [ Ixx -Ixy -Ixz ]
1127     J =| -Iyx Iyy -Iyz |
1128     [ -Izx -Iyz Izz ]
1129 chuckv 2252 */
1130    
1131     void SimInfo::getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum){
1132    
1133    
1134     double xx = 0.0;
1135     double yy = 0.0;
1136     double zz = 0.0;
1137     double xy = 0.0;
1138     double xz = 0.0;
1139     double yz = 0.0;
1140     Vector3d com(0.0);
1141     Vector3d comVel(0.0);
1142    
1143     getComAll(com, comVel);
1144    
1145     SimInfo::MoleculeIterator i;
1146     Molecule* mol;
1147    
1148     Vector3d thisq(0.0);
1149     Vector3d thisv(0.0);
1150    
1151     double thisMass = 0.0;
1152    
1153    
1154    
1155    
1156     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1157    
1158     thisq = mol->getCom()-com;
1159     thisv = mol->getComVel()-comVel;
1160     thisMass = mol->getMass();
1161     // Compute moment of intertia coefficients.
1162     xx += thisq[0]*thisq[0]*thisMass;
1163     yy += thisq[1]*thisq[1]*thisMass;
1164     zz += thisq[2]*thisq[2]*thisMass;
1165    
1166     // compute products of intertia
1167     xy += thisq[0]*thisq[1]*thisMass;
1168     xz += thisq[0]*thisq[2]*thisMass;
1169     yz += thisq[1]*thisq[2]*thisMass;
1170    
1171     angularMomentum += cross( thisq, thisv ) * thisMass;
1172    
1173     }
1174    
1175    
1176     inertiaTensor(0,0) = yy + zz;
1177     inertiaTensor(0,1) = -xy;
1178     inertiaTensor(0,2) = -xz;
1179     inertiaTensor(1,0) = -xy;
1180 chuckv 2256 inertiaTensor(1,1) = xx + zz;
1181 chuckv 2252 inertiaTensor(1,2) = -yz;
1182     inertiaTensor(2,0) = -xz;
1183     inertiaTensor(2,1) = -yz;
1184     inertiaTensor(2,2) = xx + yy;
1185    
1186     #ifdef IS_MPI
1187     Mat3x3d tmpI(inertiaTensor);
1188     Vector3d tmpAngMom;
1189     MPI_Allreduce(tmpI.getArrayPointer(), inertiaTensor.getArrayPointer(),9,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1190     MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1191     #endif
1192    
1193     return;
1194     }
1195    
1196     //Returns the angular momentum of the system
1197     Vector3d SimInfo::getAngularMomentum(){
1198    
1199     Vector3d com(0.0);
1200     Vector3d comVel(0.0);
1201     Vector3d angularMomentum(0.0);
1202    
1203     getComAll(com,comVel);
1204    
1205     SimInfo::MoleculeIterator i;
1206     Molecule* mol;
1207    
1208 chuckv 2256 Vector3d thisr(0.0);
1209     Vector3d thisp(0.0);
1210 chuckv 2252
1211 chuckv 2256 double thisMass;
1212 chuckv 2252
1213     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1214 chuckv 2256 thisMass = mol->getMass();
1215     thisr = mol->getCom()-com;
1216     thisp = (mol->getComVel()-comVel)*thisMass;
1217 chuckv 2252
1218 chuckv 2256 angularMomentum += cross( thisr, thisp );
1219    
1220 chuckv 2252 }
1221    
1222     #ifdef IS_MPI
1223     Vector3d tmpAngMom;
1224     MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1225     #endif
1226    
1227     return angularMomentum;
1228     }
1229    
1230    
1231 gezelter 1930 }//end namespace oopse
1232