ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/brains/SimInfo.cpp
Revision: 2463
Committed: Mon Nov 21 22:59:21 2005 UTC (18 years, 7 months ago) by gezelter
File size: 43289 byte(s)
Log Message:
Cutoff Mixing fixes

File Contents

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