ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/brains/SimInfo.cpp
Revision: 3013
Committed: Thu Sep 21 18:25:17 2006 UTC (17 years, 11 months ago) by chrisfen
File size: 46146 byte(s)
Log Message:
fixed the half self term for wolf electrostatics and OOPSE now chooses a cutoff radius dependent alpha for damped electrostatics

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