ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/SimInfo.cpp
Revision: 1241
Committed: Fri Apr 25 15:14:47 2008 UTC (17 years ago) by gezelter
File size: 48737 byte(s)
Log Message:
A bunch of minor changes to make MPI compilation faster than
the double compilation we do now...

File Contents

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