ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/brains/SimInfo.cpp
Revision: 2917
Committed: Mon Jul 3 13:18:43 2006 UTC (18 years, 2 months ago) by chrisfen
File size: 45290 byte(s)
Log Message:
Added simulation box dipole moment accumulation for the purposes of calculating dielectric constants

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 chuckv 2533 #include "io/ForceFieldOptions.hpp"
67     #include "UseTheForce/ForceField.hpp"
68 gezelter 1490
69 gezelter 1930 #ifdef IS_MPI
70     #include "UseTheForce/mpiComponentPlan.h"
71     #include "UseTheForce/DarkSide/simParallel_interface.h"
72     #endif
73 gezelter 1490
74 gezelter 1930 namespace oopse {
75 tim 2448 std::set<int> getRigidSet(int index, std::map<int, std::set<int> >& container) {
76     std::map<int, std::set<int> >::iterator i = container.find(index);
77     std::set<int> result;
78     if (i != container.end()) {
79     result = i->second;
80     }
81 gezelter 1490
82 tim 2448 return result;
83     }
84    
85 tim 2469 SimInfo::SimInfo(ForceField* ff, Globals* simParams) :
86     forceField_(ff), simParams_(simParams),
87 gezelter 2733 ndf_(0), fdf_local(0), ndfRaw_(0), ndfTrans_(0), nZconstraint_(0),
88 gezelter 2204 nGlobalMols_(0), nGlobalAtoms_(0), nGlobalCutoffGroups_(0),
89     nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0),
90     nAtoms_(0), nBonds_(0), nBends_(0), nTorsions_(0), nRigidBodies_(0),
91     nIntegrableObjects_(0), nCutoffGroups_(0), nConstraints_(0),
92 chrisfen 2917 sman_(NULL), fortranInitialized_(false), calcBoxDipole_(false) {
93 gezelter 1490
94 gezelter 2204 MoleculeStamp* molStamp;
95     int nMolWithSameStamp;
96     int nCutoffAtoms = 0; // number of atoms belong to cutoff groups
97 chrisfen 2344 int nGroups = 0; //total cutoff groups defined in meta-data file
98 gezelter 2204 CutoffGroupStamp* cgStamp;
99     RigidBodyStamp* rbStamp;
100     int nRigidAtoms = 0;
101 tim 2469 std::vector<Component*> components = simParams->getComponents();
102    
103     for (std::vector<Component*>::iterator i = components.begin(); i !=components.end(); ++i) {
104     molStamp = (*i)->getMoleculeStamp();
105     nMolWithSameStamp = (*i)->getNMol();
106 gezelter 1930
107     addMoleculeStamp(molStamp, nMolWithSameStamp);
108 gezelter 1490
109 gezelter 1930 //calculate atoms in molecules
110     nGlobalAtoms_ += molStamp->getNAtoms() *nMolWithSameStamp;
111 gezelter 1490
112 gezelter 1930 //calculate atoms in cutoff groups
113     int nAtomsInGroups = 0;
114     int nCutoffGroupsInStamp = molStamp->getNCutoffGroups();
115    
116     for (int j=0; j < nCutoffGroupsInStamp; j++) {
117 tim 2469 cgStamp = molStamp->getCutoffGroupStamp(j);
118 gezelter 2204 nAtomsInGroups += cgStamp->getNMembers();
119 gezelter 1930 }
120 gezelter 1490
121 gezelter 1930 nGroups += nCutoffGroupsInStamp * nMolWithSameStamp;
122 chrisfen 2344
123 gezelter 1930 nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp;
124 gezelter 1490
125 gezelter 1930 //calculate atoms in rigid bodies
126     int nAtomsInRigidBodies = 0;
127 tim 1958 int nRigidBodiesInStamp = molStamp->getNRigidBodies();
128 gezelter 1930
129     for (int j=0; j < nRigidBodiesInStamp; j++) {
130 tim 2469 rbStamp = molStamp->getRigidBodyStamp(j);
131 gezelter 2204 nAtomsInRigidBodies += rbStamp->getNMembers();
132 gezelter 1930 }
133 gezelter 1490
134 gezelter 1930 nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp;
135     nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp;
136    
137 gezelter 2204 }
138 chrisfen 1636
139 chrisfen 2344 //every free atom (atom does not belong to cutoff groups) is a cutoff
140     //group therefore the total number of cutoff groups in the system is
141     //equal to the total number of atoms minus number of atoms belong to
142     //cutoff group defined in meta-data file plus the number of cutoff
143     //groups defined in meta-data file
144 gezelter 2204 nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
145 gezelter 1490
146 chrisfen 2344 //every free atom (atom does not belong to rigid bodies) is an
147     //integrable object therefore the total number of integrable objects
148     //in the system is equal to the total number of atoms minus number of
149     //atoms belong to rigid body defined in meta-data file plus the number
150     //of rigid bodies defined in meta-data file
151     nGlobalIntegrableObjects_ = nGlobalAtoms_ - nRigidAtoms
152     + nGlobalRigidBodies_;
153    
154 gezelter 2204 nGlobalMols_ = molStampIds_.size();
155 gezelter 1490
156 gezelter 1930 #ifdef IS_MPI
157 gezelter 2204 molToProcMap_.resize(nGlobalMols_);
158 gezelter 1930 #endif
159 tim 1976
160 gezelter 2204 }
161 gezelter 1490
162 gezelter 2204 SimInfo::~SimInfo() {
163 tim 2082 std::map<int, Molecule*>::iterator i;
164     for (i = molecules_.begin(); i != molecules_.end(); ++i) {
165 gezelter 2204 delete i->second;
166 tim 2082 }
167     molecules_.clear();
168 tim 2187
169 gezelter 1930 delete sman_;
170     delete simParams_;
171     delete forceField_;
172 gezelter 2204 }
173 gezelter 1490
174 gezelter 2204 int SimInfo::getNGlobalConstraints() {
175 gezelter 1930 int nGlobalConstraints;
176     #ifdef IS_MPI
177     MPI_Allreduce(&nConstraints_, &nGlobalConstraints, 1, MPI_INT, MPI_SUM,
178     MPI_COMM_WORLD);
179     #else
180     nGlobalConstraints = nConstraints_;
181     #endif
182     return nGlobalConstraints;
183 gezelter 2204 }
184 gezelter 1490
185 gezelter 2204 bool SimInfo::addMolecule(Molecule* mol) {
186 gezelter 1930 MoleculeIterator i;
187 gezelter 1490
188 gezelter 1930 i = molecules_.find(mol->getGlobalIndex());
189     if (i == molecules_.end() ) {
190 gezelter 1490
191 gezelter 2204 molecules_.insert(std::make_pair(mol->getGlobalIndex(), mol));
192 gezelter 1930
193 gezelter 2204 nAtoms_ += mol->getNAtoms();
194     nBonds_ += mol->getNBonds();
195     nBends_ += mol->getNBends();
196     nTorsions_ += mol->getNTorsions();
197     nRigidBodies_ += mol->getNRigidBodies();
198     nIntegrableObjects_ += mol->getNIntegrableObjects();
199     nCutoffGroups_ += mol->getNCutoffGroups();
200     nConstraints_ += mol->getNConstraintPairs();
201 gezelter 1490
202 gezelter 2204 addExcludePairs(mol);
203 gezelter 1930
204 gezelter 2204 return true;
205 gezelter 1930 } else {
206 gezelter 2204 return false;
207 gezelter 1930 }
208 gezelter 2204 }
209 gezelter 1490
210 gezelter 2204 bool SimInfo::removeMolecule(Molecule* mol) {
211 gezelter 1930 MoleculeIterator i;
212     i = molecules_.find(mol->getGlobalIndex());
213 gezelter 1490
214 gezelter 1930 if (i != molecules_.end() ) {
215 gezelter 1490
216 gezelter 2204 assert(mol == i->second);
217 gezelter 1930
218 gezelter 2204 nAtoms_ -= mol->getNAtoms();
219     nBonds_ -= mol->getNBonds();
220     nBends_ -= mol->getNBends();
221     nTorsions_ -= mol->getNTorsions();
222     nRigidBodies_ -= mol->getNRigidBodies();
223     nIntegrableObjects_ -= mol->getNIntegrableObjects();
224     nCutoffGroups_ -= mol->getNCutoffGroups();
225     nConstraints_ -= mol->getNConstraintPairs();
226 gezelter 1490
227 gezelter 2204 removeExcludePairs(mol);
228     molecules_.erase(mol->getGlobalIndex());
229 gezelter 1490
230 gezelter 2204 delete mol;
231 gezelter 1930
232 gezelter 2204 return true;
233 gezelter 1930 } else {
234 gezelter 2204 return false;
235 gezelter 1930 }
236    
237    
238 gezelter 2204 }
239 gezelter 1930
240    
241 gezelter 2204 Molecule* SimInfo::beginMolecule(MoleculeIterator& i) {
242 gezelter 1930 i = molecules_.begin();
243     return i == molecules_.end() ? NULL : i->second;
244 gezelter 2204 }
245 gezelter 1930
246 gezelter 2204 Molecule* SimInfo::nextMolecule(MoleculeIterator& i) {
247 gezelter 1930 ++i;
248     return i == molecules_.end() ? NULL : i->second;
249 gezelter 2204 }
250 gezelter 1490
251    
252 gezelter 2204 void SimInfo::calcNdf() {
253 gezelter 1930 int ndf_local;
254     MoleculeIterator i;
255     std::vector<StuntDouble*>::iterator j;
256     Molecule* mol;
257     StuntDouble* integrableObject;
258 gezelter 1490
259 gezelter 1930 ndf_local = 0;
260    
261     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
262 gezelter 2204 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
263     integrableObject = mol->nextIntegrableObject(j)) {
264 gezelter 1490
265 gezelter 2204 ndf_local += 3;
266 gezelter 1490
267 gezelter 2204 if (integrableObject->isDirectional()) {
268     if (integrableObject->isLinear()) {
269     ndf_local += 2;
270     } else {
271     ndf_local += 3;
272     }
273     }
274 gezelter 1930
275 tim 2469 }
276     }
277 gezelter 1930
278     // n_constraints is local, so subtract them on each processor
279     ndf_local -= nConstraints_;
280    
281     #ifdef IS_MPI
282     MPI_Allreduce(&ndf_local,&ndf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
283     #else
284     ndf_ = ndf_local;
285     #endif
286    
287     // nZconstraints_ is global, as are the 3 COM translations for the
288     // entire system:
289     ndf_ = ndf_ - 3 - nZconstraint_;
290    
291 gezelter 2204 }
292 gezelter 1490
293 gezelter 2733 int SimInfo::getFdf() {
294     #ifdef IS_MPI
295     MPI_Allreduce(&fdf_local,&fdf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
296     #else
297     fdf_ = fdf_local;
298     #endif
299     return fdf_;
300     }
301    
302 gezelter 2204 void SimInfo::calcNdfRaw() {
303 gezelter 1930 int ndfRaw_local;
304 gezelter 1490
305 gezelter 1930 MoleculeIterator i;
306     std::vector<StuntDouble*>::iterator j;
307     Molecule* mol;
308     StuntDouble* integrableObject;
309    
310     // Raw degrees of freedom that we have to set
311     ndfRaw_local = 0;
312    
313     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
314 gezelter 2204 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
315     integrableObject = mol->nextIntegrableObject(j)) {
316 gezelter 1930
317 gezelter 2204 ndfRaw_local += 3;
318 gezelter 1930
319 gezelter 2204 if (integrableObject->isDirectional()) {
320     if (integrableObject->isLinear()) {
321     ndfRaw_local += 2;
322     } else {
323     ndfRaw_local += 3;
324     }
325     }
326 gezelter 1930
327 gezelter 2204 }
328 gezelter 1930 }
329    
330     #ifdef IS_MPI
331     MPI_Allreduce(&ndfRaw_local,&ndfRaw_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
332     #else
333     ndfRaw_ = ndfRaw_local;
334     #endif
335 gezelter 2204 }
336 gezelter 1490
337 gezelter 2204 void SimInfo::calcNdfTrans() {
338 gezelter 1930 int ndfTrans_local;
339 gezelter 1490
340 gezelter 1930 ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_;
341 gezelter 1490
342    
343 gezelter 1930 #ifdef IS_MPI
344     MPI_Allreduce(&ndfTrans_local,&ndfTrans_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
345     #else
346     ndfTrans_ = ndfTrans_local;
347     #endif
348 gezelter 1490
349 gezelter 1930 ndfTrans_ = ndfTrans_ - 3 - nZconstraint_;
350    
351 gezelter 2204 }
352 gezelter 1490
353 gezelter 2204 void SimInfo::addExcludePairs(Molecule* mol) {
354 gezelter 1930 std::vector<Bond*>::iterator bondIter;
355     std::vector<Bend*>::iterator bendIter;
356     std::vector<Torsion*>::iterator torsionIter;
357     Bond* bond;
358     Bend* bend;
359     Torsion* torsion;
360     int a;
361     int b;
362     int c;
363     int d;
364 tim 2448
365     std::map<int, std::set<int> > atomGroups;
366    
367     Molecule::RigidBodyIterator rbIter;
368     RigidBody* rb;
369     Molecule::IntegrableObjectIterator ii;
370     StuntDouble* integrableObject;
371 gezelter 1930
372 tim 2448 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
373     integrableObject = mol->nextIntegrableObject(ii)) {
374    
375     if (integrableObject->isRigidBody()) {
376     rb = static_cast<RigidBody*>(integrableObject);
377     std::vector<Atom*> atoms = rb->getAtoms();
378     std::set<int> rigidAtoms;
379     for (int i = 0; i < atoms.size(); ++i) {
380     rigidAtoms.insert(atoms[i]->getGlobalIndex());
381     }
382     for (int i = 0; i < atoms.size(); ++i) {
383     atomGroups.insert(std::map<int, std::set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
384     }
385     } else {
386     std::set<int> oneAtomSet;
387     oneAtomSet.insert(integrableObject->getGlobalIndex());
388     atomGroups.insert(std::map<int, std::set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));
389     }
390     }
391    
392    
393    
394 gezelter 1930 for (bond= mol->beginBond(bondIter); bond != NULL; bond = mol->nextBond(bondIter)) {
395 gezelter 2204 a = bond->getAtomA()->getGlobalIndex();
396     b = bond->getAtomB()->getGlobalIndex();
397     exclude_.addPair(a, b);
398 gezelter 1930 }
399 gezelter 1490
400 gezelter 1930 for (bend= mol->beginBend(bendIter); bend != NULL; bend = mol->nextBend(bendIter)) {
401 gezelter 2204 a = bend->getAtomA()->getGlobalIndex();
402     b = bend->getAtomB()->getGlobalIndex();
403     c = bend->getAtomC()->getGlobalIndex();
404 tim 2448 std::set<int> rigidSetA = getRigidSet(a, atomGroups);
405     std::set<int> rigidSetB = getRigidSet(b, atomGroups);
406     std::set<int> rigidSetC = getRigidSet(c, atomGroups);
407 gezelter 1490
408 tim 2448 exclude_.addPairs(rigidSetA, rigidSetB);
409     exclude_.addPairs(rigidSetA, rigidSetC);
410     exclude_.addPairs(rigidSetB, rigidSetC);
411    
412     //exclude_.addPair(a, b);
413     //exclude_.addPair(a, c);
414     //exclude_.addPair(b, c);
415 gezelter 1930 }
416 gezelter 1490
417 gezelter 1930 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL; torsion = mol->nextTorsion(torsionIter)) {
418 gezelter 2204 a = torsion->getAtomA()->getGlobalIndex();
419     b = torsion->getAtomB()->getGlobalIndex();
420     c = torsion->getAtomC()->getGlobalIndex();
421     d = torsion->getAtomD()->getGlobalIndex();
422 tim 2448 std::set<int> rigidSetA = getRigidSet(a, atomGroups);
423     std::set<int> rigidSetB = getRigidSet(b, atomGroups);
424     std::set<int> rigidSetC = getRigidSet(c, atomGroups);
425     std::set<int> rigidSetD = getRigidSet(d, atomGroups);
426 gezelter 1490
427 tim 2448 exclude_.addPairs(rigidSetA, rigidSetB);
428     exclude_.addPairs(rigidSetA, rigidSetC);
429     exclude_.addPairs(rigidSetA, rigidSetD);
430     exclude_.addPairs(rigidSetB, rigidSetC);
431     exclude_.addPairs(rigidSetB, rigidSetD);
432     exclude_.addPairs(rigidSetC, rigidSetD);
433    
434     /*
435     exclude_.addPairs(rigidSetA.begin(), rigidSetA.end(), rigidSetB.begin(), rigidSetB.end());
436     exclude_.addPairs(rigidSetA.begin(), rigidSetA.end(), rigidSetC.begin(), rigidSetC.end());
437     exclude_.addPairs(rigidSetA.begin(), rigidSetA.end(), rigidSetD.begin(), rigidSetD.end());
438     exclude_.addPairs(rigidSetB.begin(), rigidSetB.end(), rigidSetC.begin(), rigidSetC.end());
439     exclude_.addPairs(rigidSetB.begin(), rigidSetB.end(), rigidSetD.begin(), rigidSetD.end());
440     exclude_.addPairs(rigidSetC.begin(), rigidSetC.end(), rigidSetD.begin(), rigidSetD.end());
441    
442    
443 gezelter 2204 exclude_.addPair(a, b);
444     exclude_.addPair(a, c);
445     exclude_.addPair(a, d);
446     exclude_.addPair(b, c);
447     exclude_.addPair(b, d);
448     exclude_.addPair(c, d);
449 tim 2448 */
450 gezelter 1490 }
451    
452 tim 2114 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
453 gezelter 2204 std::vector<Atom*> atoms = rb->getAtoms();
454     for (int i = 0; i < atoms.size() -1 ; ++i) {
455     for (int j = i + 1; j < atoms.size(); ++j) {
456     a = atoms[i]->getGlobalIndex();
457     b = atoms[j]->getGlobalIndex();
458     exclude_.addPair(a, b);
459     }
460     }
461 tim 2114 }
462    
463 gezelter 2204 }
464 gezelter 1930
465 gezelter 2204 void SimInfo::removeExcludePairs(Molecule* mol) {
466 gezelter 1930 std::vector<Bond*>::iterator bondIter;
467     std::vector<Bend*>::iterator bendIter;
468     std::vector<Torsion*>::iterator torsionIter;
469     Bond* bond;
470     Bend* bend;
471     Torsion* torsion;
472     int a;
473     int b;
474     int c;
475     int d;
476 tim 2448
477     std::map<int, std::set<int> > atomGroups;
478    
479     Molecule::RigidBodyIterator rbIter;
480     RigidBody* rb;
481     Molecule::IntegrableObjectIterator ii;
482     StuntDouble* integrableObject;
483 gezelter 1930
484 tim 2448 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
485     integrableObject = mol->nextIntegrableObject(ii)) {
486    
487     if (integrableObject->isRigidBody()) {
488     rb = static_cast<RigidBody*>(integrableObject);
489     std::vector<Atom*> atoms = rb->getAtoms();
490     std::set<int> rigidAtoms;
491     for (int i = 0; i < atoms.size(); ++i) {
492     rigidAtoms.insert(atoms[i]->getGlobalIndex());
493     }
494     for (int i = 0; i < atoms.size(); ++i) {
495     atomGroups.insert(std::map<int, std::set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
496     }
497     } else {
498     std::set<int> oneAtomSet;
499     oneAtomSet.insert(integrableObject->getGlobalIndex());
500     atomGroups.insert(std::map<int, std::set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));
501     }
502     }
503    
504    
505 gezelter 1930 for (bond= mol->beginBond(bondIter); bond != NULL; bond = mol->nextBond(bondIter)) {
506 gezelter 2204 a = bond->getAtomA()->getGlobalIndex();
507     b = bond->getAtomB()->getGlobalIndex();
508     exclude_.removePair(a, b);
509 gezelter 1490 }
510 gezelter 1930
511     for (bend= mol->beginBend(bendIter); bend != NULL; bend = mol->nextBend(bendIter)) {
512 gezelter 2204 a = bend->getAtomA()->getGlobalIndex();
513     b = bend->getAtomB()->getGlobalIndex();
514     c = bend->getAtomC()->getGlobalIndex();
515 gezelter 1930
516 tim 2448 std::set<int> rigidSetA = getRigidSet(a, atomGroups);
517     std::set<int> rigidSetB = getRigidSet(b, atomGroups);
518     std::set<int> rigidSetC = getRigidSet(c, atomGroups);
519    
520     exclude_.removePairs(rigidSetA, rigidSetB);
521     exclude_.removePairs(rigidSetA, rigidSetC);
522     exclude_.removePairs(rigidSetB, rigidSetC);
523    
524     //exclude_.removePair(a, b);
525     //exclude_.removePair(a, c);
526     //exclude_.removePair(b, c);
527 gezelter 1490 }
528 gezelter 1930
529     for (torsion= mol->beginTorsion(torsionIter); torsion != NULL; torsion = mol->nextTorsion(torsionIter)) {
530 gezelter 2204 a = torsion->getAtomA()->getGlobalIndex();
531     b = torsion->getAtomB()->getGlobalIndex();
532     c = torsion->getAtomC()->getGlobalIndex();
533     d = torsion->getAtomD()->getGlobalIndex();
534 gezelter 1930
535 tim 2448 std::set<int> rigidSetA = getRigidSet(a, atomGroups);
536     std::set<int> rigidSetB = getRigidSet(b, atomGroups);
537     std::set<int> rigidSetC = getRigidSet(c, atomGroups);
538     std::set<int> rigidSetD = getRigidSet(d, atomGroups);
539    
540     exclude_.removePairs(rigidSetA, rigidSetB);
541     exclude_.removePairs(rigidSetA, rigidSetC);
542     exclude_.removePairs(rigidSetA, rigidSetD);
543     exclude_.removePairs(rigidSetB, rigidSetC);
544     exclude_.removePairs(rigidSetB, rigidSetD);
545     exclude_.removePairs(rigidSetC, rigidSetD);
546    
547     /*
548     exclude_.removePairs(rigidSetA.begin(), rigidSetA.end(), rigidSetB.begin(), rigidSetB.end());
549     exclude_.removePairs(rigidSetA.begin(), rigidSetA.end(), rigidSetC.begin(), rigidSetC.end());
550     exclude_.removePairs(rigidSetA.begin(), rigidSetA.end(), rigidSetD.begin(), rigidSetD.end());
551     exclude_.removePairs(rigidSetB.begin(), rigidSetB.end(), rigidSetC.begin(), rigidSetC.end());
552     exclude_.removePairs(rigidSetB.begin(), rigidSetB.end(), rigidSetD.begin(), rigidSetD.end());
553     exclude_.removePairs(rigidSetC.begin(), rigidSetC.end(), rigidSetD.begin(), rigidSetD.end());
554    
555    
556 gezelter 2204 exclude_.removePair(a, b);
557     exclude_.removePair(a, c);
558     exclude_.removePair(a, d);
559     exclude_.removePair(b, c);
560     exclude_.removePair(b, d);
561     exclude_.removePair(c, d);
562 tim 2448 */
563 gezelter 1930 }
564    
565 tim 2114 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
566 gezelter 2204 std::vector<Atom*> atoms = rb->getAtoms();
567     for (int i = 0; i < atoms.size() -1 ; ++i) {
568     for (int j = i + 1; j < atoms.size(); ++j) {
569     a = atoms[i]->getGlobalIndex();
570     b = atoms[j]->getGlobalIndex();
571     exclude_.removePair(a, b);
572     }
573     }
574 tim 2114 }
575    
576 gezelter 2204 }
577 gezelter 1490
578    
579 gezelter 2204 void SimInfo::addMoleculeStamp(MoleculeStamp* molStamp, int nmol) {
580 gezelter 1930 int curStampId;
581 gezelter 1490
582 gezelter 1930 //index from 0
583     curStampId = moleculeStamps_.size();
584 gezelter 1490
585 gezelter 1930 moleculeStamps_.push_back(molStamp);
586     molStampIds_.insert(molStampIds_.end(), nmol, curStampId);
587 gezelter 2204 }
588 gezelter 1490
589 gezelter 2204 void SimInfo::update() {
590 gezelter 1490
591 gezelter 1930 setupSimType();
592 gezelter 1490
593 gezelter 1930 #ifdef IS_MPI
594     setupFortranParallel();
595     #endif
596 gezelter 1490
597 gezelter 1930 setupFortranSim();
598 gezelter 1490
599 gezelter 1930 //setup fortran force field
600     /** @deprecate */
601     int isError = 0;
602 chrisfen 2297
603 chrisfen 2302 setupElectrostaticSummationMethod( isError );
604 chrisfen 2425 setupSwitchingFunction();
605 chrisfen 2917 setupAccumulateBoxDipole();
606 chrisfen 2297
607 gezelter 1930 if(isError){
608 gezelter 2204 sprintf( painCave.errMsg,
609     "ForceField error: There was an error initializing the forceField in fortran.\n" );
610     painCave.isFatal = 1;
611     simError();
612 gezelter 1930 }
613 gezelter 1490
614 gezelter 1930
615     setupCutoff();
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 2297
1077     errorOut = isError;
1078 chrisfen 2309 alphaVal = simParams_->getDampingAlpha();
1079     dielectric = simParams_->getDielectric();
1080 chrisfen 2297
1081 chrisfen 2302 if (simParams_->haveElectrostaticSummationMethod()) {
1082 chrisfen 2303 std::string myMethod = simParams_->getElectrostaticSummationMethod();
1083 tim 2364 toUpper(myMethod);
1084 chrisfen 2302 if (myMethod == "NONE") {
1085     esm = NONE;
1086 chrisfen 2297 } else {
1087 chrisfen 2408 if (myMethod == "SWITCHING_FUNCTION") {
1088     esm = SWITCHING_FUNCTION;
1089 chrisfen 2297 } else {
1090 chrisfen 2408 if (myMethod == "SHIFTED_POTENTIAL") {
1091     esm = SHIFTED_POTENTIAL;
1092     } else {
1093     if (myMethod == "SHIFTED_FORCE") {
1094     esm = SHIFTED_FORCE;
1095 chrisfen 2297 } else {
1096 chrisfen 2408 if (myMethod == "REACTION_FIELD") {
1097     esm = REACTION_FIELD;
1098     } else {
1099     // throw error
1100     sprintf( painCave.errMsg,
1101 gezelter 2463 "SimInfo error: Unknown electrostaticSummationMethod.\n"
1102     "\t(Input file specified %s .)\n"
1103     "\telectrostaticSummationMethod must be one of: \"none\",\n"
1104     "\t\"shifted_potential\", \"shifted_force\", or \n"
1105     "\t\"reaction_field\".\n", myMethod.c_str() );
1106 chrisfen 2408 painCave.isFatal = 1;
1107     simError();
1108     }
1109     }
1110     }
1111 chrisfen 2297 }
1112     }
1113     }
1114 chrisfen 2408
1115 chrisfen 2415 if (simParams_->haveElectrostaticScreeningMethod()) {
1116     std::string myScreen = simParams_->getElectrostaticScreeningMethod();
1117 chrisfen 2408 toUpper(myScreen);
1118     if (myScreen == "UNDAMPED") {
1119     sm = UNDAMPED;
1120     } else {
1121     if (myScreen == "DAMPED") {
1122     sm = DAMPED;
1123     if (!simParams_->haveDampingAlpha()) {
1124     //throw error
1125     sprintf( painCave.errMsg,
1126 gezelter 2463 "SimInfo warning: dampingAlpha was not specified in the input file.\n"
1127     "\tA default value of %f (1/ang) will be used.\n", alphaVal);
1128 chrisfen 2408 painCave.isFatal = 0;
1129     simError();
1130     }
1131 chrisfen 2415 } else {
1132     // throw error
1133     sprintf( painCave.errMsg,
1134 gezelter 2463 "SimInfo error: Unknown electrostaticScreeningMethod.\n"
1135     "\t(Input file specified %s .)\n"
1136     "\telectrostaticScreeningMethod must be one of: \"undamped\"\n"
1137     "or \"damped\".\n", myScreen.c_str() );
1138 chrisfen 2415 painCave.isFatal = 1;
1139     simError();
1140 chrisfen 2408 }
1141     }
1142     }
1143 chrisfen 2415
1144 chrisfen 2309 // let's pass some summation method variables to fortran
1145 chrisfen 2552 setElectrostaticSummationMethod( &esm );
1146 gezelter 2508 setFortranElectrostaticMethod( &esm );
1147 chrisfen 2408 setScreeningMethod( &sm );
1148     setDampingAlpha( &alphaVal );
1149 chrisfen 2309 setReactionFieldDielectric( &dielectric );
1150 gezelter 2463 initFortranFF( &errorOut );
1151 chrisfen 2297 }
1152    
1153 chrisfen 2425 void SimInfo::setupSwitchingFunction() {
1154     int ft = CUBIC;
1155    
1156     if (simParams_->haveSwitchingFunctionType()) {
1157     std::string funcType = simParams_->getSwitchingFunctionType();
1158     toUpper(funcType);
1159     if (funcType == "CUBIC") {
1160     ft = CUBIC;
1161     } else {
1162     if (funcType == "FIFTH_ORDER_POLYNOMIAL") {
1163     ft = FIFTH_ORDER_POLY;
1164     } else {
1165     // throw error
1166     sprintf( painCave.errMsg,
1167     "SimInfo error: Unknown switchingFunctionType. (Input file specified %s .)\n\tswitchingFunctionType must be one of: \"cubic\" or \"fifth_order_polynomial\".", funcType.c_str() );
1168     painCave.isFatal = 1;
1169     simError();
1170     }
1171     }
1172     }
1173    
1174     // send switching function notification to switcheroo
1175     setFunctionType(&ft);
1176    
1177     }
1178    
1179 chrisfen 2917 void SimInfo::setupAccumulateBoxDipole() {
1180    
1181     // we only call setAccumulateBoxDipole if the accumulateBoxDipole parameter is true
1182     if ( simParams_->haveAccumulateBoxDipole() )
1183     if ( simParams_->getAccumulateBoxDipole() ) {
1184     setAccumulateBoxDipole();
1185     calcBoxDipole_ = true;
1186     }
1187    
1188     }
1189    
1190 gezelter 2204 void SimInfo::addProperty(GenericData* genData) {
1191 gezelter 1930 properties_.addProperty(genData);
1192 gezelter 2204 }
1193 gezelter 1490
1194 gezelter 2204 void SimInfo::removeProperty(const std::string& propName) {
1195 gezelter 1930 properties_.removeProperty(propName);
1196 gezelter 2204 }
1197 gezelter 1490
1198 gezelter 2204 void SimInfo::clearProperties() {
1199 gezelter 1930 properties_.clearProperties();
1200 gezelter 2204 }
1201 gezelter 1490
1202 gezelter 2204 std::vector<std::string> SimInfo::getPropertyNames() {
1203 gezelter 1930 return properties_.getPropertyNames();
1204 gezelter 2204 }
1205 gezelter 1930
1206 gezelter 2204 std::vector<GenericData*> SimInfo::getProperties() {
1207 gezelter 1930 return properties_.getProperties();
1208 gezelter 2204 }
1209 gezelter 1490
1210 gezelter 2204 GenericData* SimInfo::getPropertyByName(const std::string& propName) {
1211 gezelter 1930 return properties_.getPropertyByName(propName);
1212 gezelter 2204 }
1213 gezelter 1490
1214 gezelter 2204 void SimInfo::setSnapshotManager(SnapshotManager* sman) {
1215 tim 2116 if (sman_ == sman) {
1216 gezelter 2204 return;
1217 tim 2116 }
1218     delete sman_;
1219 gezelter 1930 sman_ = sman;
1220 gezelter 1490
1221 gezelter 1930 Molecule* mol;
1222     RigidBody* rb;
1223     Atom* atom;
1224     SimInfo::MoleculeIterator mi;
1225     Molecule::RigidBodyIterator rbIter;
1226     Molecule::AtomIterator atomIter;;
1227    
1228     for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
1229    
1230 gezelter 2204 for (atom = mol->beginAtom(atomIter); atom != NULL; atom = mol->nextAtom(atomIter)) {
1231     atom->setSnapshotManager(sman_);
1232     }
1233 gezelter 1930
1234 gezelter 2204 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
1235     rb->setSnapshotManager(sman_);
1236     }
1237 gezelter 1930 }
1238 gezelter 1490
1239 gezelter 2204 }
1240 gezelter 1490
1241 gezelter 2204 Vector3d SimInfo::getComVel(){
1242 gezelter 1930 SimInfo::MoleculeIterator i;
1243     Molecule* mol;
1244 gezelter 1490
1245 gezelter 1930 Vector3d comVel(0.0);
1246 tim 2759 RealType totalMass = 0.0;
1247 gezelter 1490
1248 gezelter 1930
1249     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1250 tim 2759 RealType mass = mol->getMass();
1251 gezelter 2204 totalMass += mass;
1252     comVel += mass * mol->getComVel();
1253 gezelter 1930 }
1254 gezelter 1490
1255 gezelter 1930 #ifdef IS_MPI
1256 tim 2759 RealType tmpMass = totalMass;
1257 gezelter 1930 Vector3d tmpComVel(comVel);
1258 tim 2759 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1259     MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1260 gezelter 1930 #endif
1261    
1262     comVel /= totalMass;
1263    
1264     return comVel;
1265 gezelter 2204 }
1266 gezelter 1490
1267 gezelter 2204 Vector3d SimInfo::getCom(){
1268 gezelter 1930 SimInfo::MoleculeIterator i;
1269     Molecule* mol;
1270 gezelter 1490
1271 gezelter 1930 Vector3d com(0.0);
1272 tim 2759 RealType totalMass = 0.0;
1273 gezelter 1930
1274     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1275 tim 2759 RealType mass = mol->getMass();
1276 gezelter 2204 totalMass += mass;
1277     com += mass * mol->getCom();
1278 gezelter 1930 }
1279 gezelter 1490
1280     #ifdef IS_MPI
1281 tim 2759 RealType tmpMass = totalMass;
1282 gezelter 1930 Vector3d tmpCom(com);
1283 tim 2759 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1284     MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1285 gezelter 1490 #endif
1286    
1287 gezelter 1930 com /= totalMass;
1288 gezelter 1490
1289 gezelter 1930 return com;
1290 gezelter 1490
1291 gezelter 2204 }
1292 gezelter 1930
1293 gezelter 2204 std::ostream& operator <<(std::ostream& o, SimInfo& info) {
1294 gezelter 1930
1295     return o;
1296 gezelter 2204 }
1297 chuckv 2252
1298    
1299     /*
1300     Returns center of mass and center of mass velocity in one function call.
1301     */
1302    
1303     void SimInfo::getComAll(Vector3d &com, Vector3d &comVel){
1304     SimInfo::MoleculeIterator i;
1305     Molecule* mol;
1306    
1307    
1308 tim 2759 RealType totalMass = 0.0;
1309 chuckv 2252
1310 gezelter 1930
1311 chuckv 2252 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1312 tim 2759 RealType mass = mol->getMass();
1313 chuckv 2252 totalMass += mass;
1314     com += mass * mol->getCom();
1315     comVel += mass * mol->getComVel();
1316     }
1317    
1318     #ifdef IS_MPI
1319 tim 2759 RealType tmpMass = totalMass;
1320 chuckv 2252 Vector3d tmpCom(com);
1321     Vector3d tmpComVel(comVel);
1322 tim 2759 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1323     MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1324     MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1325 chuckv 2252 #endif
1326    
1327     com /= totalMass;
1328     comVel /= totalMass;
1329     }
1330    
1331     /*
1332     Return intertia tensor for entire system and angular momentum Vector.
1333 chuckv 2256
1334    
1335     [ Ixx -Ixy -Ixz ]
1336     J =| -Iyx Iyy -Iyz |
1337     [ -Izx -Iyz Izz ]
1338 chuckv 2252 */
1339    
1340     void SimInfo::getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum){
1341    
1342    
1343 tim 2759 RealType xx = 0.0;
1344     RealType yy = 0.0;
1345     RealType zz = 0.0;
1346     RealType xy = 0.0;
1347     RealType xz = 0.0;
1348     RealType yz = 0.0;
1349 chuckv 2252 Vector3d com(0.0);
1350     Vector3d comVel(0.0);
1351    
1352     getComAll(com, comVel);
1353    
1354     SimInfo::MoleculeIterator i;
1355     Molecule* mol;
1356    
1357     Vector3d thisq(0.0);
1358     Vector3d thisv(0.0);
1359    
1360 tim 2759 RealType thisMass = 0.0;
1361 chuckv 2252
1362    
1363    
1364    
1365     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1366    
1367     thisq = mol->getCom()-com;
1368     thisv = mol->getComVel()-comVel;
1369     thisMass = mol->getMass();
1370     // Compute moment of intertia coefficients.
1371     xx += thisq[0]*thisq[0]*thisMass;
1372     yy += thisq[1]*thisq[1]*thisMass;
1373     zz += thisq[2]*thisq[2]*thisMass;
1374    
1375     // compute products of intertia
1376     xy += thisq[0]*thisq[1]*thisMass;
1377     xz += thisq[0]*thisq[2]*thisMass;
1378     yz += thisq[1]*thisq[2]*thisMass;
1379    
1380     angularMomentum += cross( thisq, thisv ) * thisMass;
1381    
1382     }
1383    
1384    
1385     inertiaTensor(0,0) = yy + zz;
1386     inertiaTensor(0,1) = -xy;
1387     inertiaTensor(0,2) = -xz;
1388     inertiaTensor(1,0) = -xy;
1389 chuckv 2256 inertiaTensor(1,1) = xx + zz;
1390 chuckv 2252 inertiaTensor(1,2) = -yz;
1391     inertiaTensor(2,0) = -xz;
1392     inertiaTensor(2,1) = -yz;
1393     inertiaTensor(2,2) = xx + yy;
1394    
1395     #ifdef IS_MPI
1396     Mat3x3d tmpI(inertiaTensor);
1397     Vector3d tmpAngMom;
1398 tim 2759 MPI_Allreduce(tmpI.getArrayPointer(), inertiaTensor.getArrayPointer(),9,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1399     MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1400 chuckv 2252 #endif
1401    
1402     return;
1403     }
1404    
1405     //Returns the angular momentum of the system
1406     Vector3d SimInfo::getAngularMomentum(){
1407    
1408     Vector3d com(0.0);
1409     Vector3d comVel(0.0);
1410     Vector3d angularMomentum(0.0);
1411    
1412     getComAll(com,comVel);
1413    
1414     SimInfo::MoleculeIterator i;
1415     Molecule* mol;
1416    
1417 chuckv 2256 Vector3d thisr(0.0);
1418     Vector3d thisp(0.0);
1419 chuckv 2252
1420 tim 2759 RealType thisMass;
1421 chuckv 2252
1422     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1423 chuckv 2256 thisMass = mol->getMass();
1424     thisr = mol->getCom()-com;
1425     thisp = (mol->getComVel()-comVel)*thisMass;
1426 chuckv 2252
1427 chuckv 2256 angularMomentum += cross( thisr, thisp );
1428    
1429 chuckv 2252 }
1430    
1431     #ifdef IS_MPI
1432     Vector3d tmpAngMom;
1433 tim 2759 MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1434 chuckv 2252 #endif
1435    
1436     return angularMomentum;
1437     }
1438    
1439    
1440 gezelter 1930 }//end namespace oopse
1441