ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/brains/SimInfo.cpp
Revision: 2433
Committed: Tue Nov 15 16:05:38 2005 UTC (18 years, 7 months ago) by chuckv
File size: 37827 byte(s)
Log Message:
Sutton-Chen added to SimInfo

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