ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/brains/SimInfo.cpp
Revision: 2469
Committed: Fri Dec 2 15:38:03 2005 UTC (18 years, 7 months ago) by tim
File size: 43176 byte(s)
Log Message:
End of the Link --> List
Return of the Oject-Oriented
replace yacc/lex parser with antlr parser

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