ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.cpp
Revision: 1534
Committed: Wed Dec 29 21:53:28 2010 UTC (14 years, 4 months ago) by gezelter
File size: 44577 byte(s)
Log Message:
getting closer to compiling

File Contents

# User Rev Content
1 gezelter 507 /*
2 gezelter 246 * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3     *
4     * The University of Notre Dame grants you ("Licensee") a
5     * non-exclusive, royalty free, license to use, modify and
6     * redistribute this software in source and binary code form, provided
7     * that the following conditions are met:
8     *
9 gezelter 1390 * 1. Redistributions of source code must retain the above copyright
10 gezelter 246 * notice, this list of conditions and the following disclaimer.
11     *
12 gezelter 1390 * 2. Redistributions in binary form must reproduce the above copyright
13 gezelter 246 * notice, this list of conditions and the following disclaimer in the
14     * documentation and/or other materials provided with the
15     * distribution.
16     *
17     * This software is provided "AS IS," without a warranty of any
18     * kind. All express or implied conditions, representations and
19     * warranties, including any implied warranty of merchantability,
20     * fitness for a particular purpose or non-infringement, are hereby
21     * excluded. The University of Notre Dame and its licensors shall not
22     * be liable for any damages suffered by licensee as a result of
23     * using, modifying or distributing the software or its
24     * derivatives. In no event will the University of Notre Dame or its
25     * licensors be liable for any lost revenue, profit or data, or for
26     * direct, indirect, special, consequential, incidental or punitive
27     * damages, however caused and regardless of the theory of liability,
28     * arising out of the use of or inability to use software, even if the
29     * University of Notre Dame has been advised of the possibility of
30     * such damages.
31 gezelter 1390 *
32     * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33     * research, please cite the appropriate papers when you publish your
34     * work. Good starting points are:
35     *
36     * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37     * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38     * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).
39     * [4] Vardeman & Gezelter, in progress (2009).
40 gezelter 246 */
41    
42     /**
43     * @file SimInfo.cpp
44     * @author tlin
45     * @date 11/02/2004
46     * @version 1.0
47     */
48 gezelter 2
49 gezelter 246 #include <algorithm>
50     #include <set>
51 tim 749 #include <map>
52 gezelter 2
53 tim 3 #include "brains/SimInfo.hpp"
54 gezelter 246 #include "math/Vector3.hpp"
55     #include "primitives/Molecule.hpp"
56 tim 1024 #include "primitives/StuntDouble.hpp"
57 gezelter 246 #include "UseTheForce/doForces_interface.h"
58 chuckv 1095 #include "UseTheForce/DarkSide/neighborLists_interface.h"
59 gezelter 246 #include "utils/MemoryUtils.hpp"
60 tim 3 #include "utils/simError.h"
61 tim 316 #include "selection/SelectionManager.hpp"
62 chuckv 834 #include "io/ForceFieldOptions.hpp"
63     #include "UseTheForce/ForceField.hpp"
64 gezelter 1530 #include "nonbonded/SwitchingFunction.hpp"
65 gezelter 2
66 gezelter 246 #ifdef IS_MPI
67     #include "UseTheForce/mpiComponentPlan.h"
68     #include "UseTheForce/DarkSide/simParallel_interface.h"
69     #endif
70 gezelter 2
71 gezelter 1528 using namespace std;
72 gezelter 1390 namespace OpenMD {
73 tim 749
74 tim 770 SimInfo::SimInfo(ForceField* ff, Globals* simParams) :
75     forceField_(ff), simParams_(simParams),
76 gezelter 945 ndf_(0), fdf_local(0), ndfRaw_(0), ndfTrans_(0), nZconstraint_(0),
77 gezelter 507 nGlobalMols_(0), nGlobalAtoms_(0), nGlobalCutoffGroups_(0),
78     nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0),
79 gezelter 1277 nAtoms_(0), nBonds_(0), nBends_(0), nTorsions_(0), nInversions_(0),
80     nRigidBodies_(0), nIntegrableObjects_(0), nCutoffGroups_(0),
81     nConstraints_(0), sman_(NULL), fortranInitialized_(false),
82 gezelter 1528 calcBoxDipole_(false), useAtomicVirial_(true) {
83    
84     MoleculeStamp* molStamp;
85     int nMolWithSameStamp;
86     int nCutoffAtoms = 0; // number of atoms belong to cutoff groups
87     int nGroups = 0; //total cutoff groups defined in meta-data file
88     CutoffGroupStamp* cgStamp;
89     RigidBodyStamp* rbStamp;
90     int nRigidAtoms = 0;
91    
92     vector<Component*> components = simParams->getComponents();
93    
94     for (vector<Component*>::iterator i = components.begin(); i !=components.end(); ++i) {
95     molStamp = (*i)->getMoleculeStamp();
96     nMolWithSameStamp = (*i)->getNMol();
97 tim 770
98 gezelter 1528 addMoleculeStamp(molStamp, nMolWithSameStamp);
99    
100     //calculate atoms in molecules
101     nGlobalAtoms_ += molStamp->getNAtoms() *nMolWithSameStamp;
102    
103     //calculate atoms in cutoff groups
104     int nAtomsInGroups = 0;
105     int nCutoffGroupsInStamp = molStamp->getNCutoffGroups();
106    
107     for (int j=0; j < nCutoffGroupsInStamp; j++) {
108     cgStamp = molStamp->getCutoffGroupStamp(j);
109     nAtomsInGroups += cgStamp->getNMembers();
110 gezelter 507 }
111 gezelter 1528
112     nGroups += nCutoffGroupsInStamp * nMolWithSameStamp;
113    
114     nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp;
115    
116     //calculate atoms in rigid bodies
117     int nAtomsInRigidBodies = 0;
118     int nRigidBodiesInStamp = molStamp->getNRigidBodies();
119    
120     for (int j=0; j < nRigidBodiesInStamp; j++) {
121     rbStamp = molStamp->getRigidBodyStamp(j);
122     nAtomsInRigidBodies += rbStamp->getNMembers();
123     }
124    
125     nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp;
126     nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp;
127    
128     }
129    
130     //every free atom (atom does not belong to cutoff groups) is a cutoff
131     //group therefore the total number of cutoff groups in the system is
132     //equal to the total number of atoms minus number of atoms belong to
133     //cutoff group defined in meta-data file plus the number of cutoff
134     //groups defined in meta-data file
135     nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
136    
137     //every free atom (atom does not belong to rigid bodies) is an
138     //integrable object therefore the total number of integrable objects
139     //in the system is equal to the total number of atoms minus number of
140     //atoms belong to rigid body defined in meta-data file plus the number
141     //of rigid bodies defined in meta-data file
142     nGlobalIntegrableObjects_ = nGlobalAtoms_ - nRigidAtoms
143     + nGlobalRigidBodies_;
144    
145     nGlobalMols_ = molStampIds_.size();
146     molToProcMap_.resize(nGlobalMols_);
147     }
148 chrisfen 645
149 gezelter 507 SimInfo::~SimInfo() {
150 gezelter 1528 map<int, Molecule*>::iterator i;
151 tim 398 for (i = molecules_.begin(); i != molecules_.end(); ++i) {
152 gezelter 507 delete i->second;
153 tim 398 }
154     molecules_.clear();
155 tim 490
156 gezelter 246 delete sman_;
157     delete simParams_;
158     delete forceField_;
159 gezelter 507 }
160 gezelter 2
161    
162 gezelter 507 bool SimInfo::addMolecule(Molecule* mol) {
163 gezelter 246 MoleculeIterator i;
164 gezelter 1528
165 gezelter 246 i = molecules_.find(mol->getGlobalIndex());
166     if (i == molecules_.end() ) {
167 gezelter 1528
168     molecules_.insert(make_pair(mol->getGlobalIndex(), mol));
169    
170 gezelter 507 nAtoms_ += mol->getNAtoms();
171     nBonds_ += mol->getNBonds();
172     nBends_ += mol->getNBends();
173     nTorsions_ += mol->getNTorsions();
174 gezelter 1277 nInversions_ += mol->getNInversions();
175 gezelter 507 nRigidBodies_ += mol->getNRigidBodies();
176     nIntegrableObjects_ += mol->getNIntegrableObjects();
177     nCutoffGroups_ += mol->getNCutoffGroups();
178     nConstraints_ += mol->getNConstraintPairs();
179 gezelter 1528
180 gezelter 1287 addInteractionPairs(mol);
181 gezelter 1528
182 gezelter 507 return true;
183 gezelter 246 } else {
184 gezelter 507 return false;
185 gezelter 246 }
186 gezelter 507 }
187 gezelter 1528
188 gezelter 507 bool SimInfo::removeMolecule(Molecule* mol) {
189 gezelter 246 MoleculeIterator i;
190     i = molecules_.find(mol->getGlobalIndex());
191 gezelter 2
192 gezelter 246 if (i != molecules_.end() ) {
193 gezelter 2
194 gezelter 507 assert(mol == i->second);
195 gezelter 246
196 gezelter 507 nAtoms_ -= mol->getNAtoms();
197     nBonds_ -= mol->getNBonds();
198     nBends_ -= mol->getNBends();
199     nTorsions_ -= mol->getNTorsions();
200 gezelter 1277 nInversions_ -= mol->getNInversions();
201 gezelter 507 nRigidBodies_ -= mol->getNRigidBodies();
202     nIntegrableObjects_ -= mol->getNIntegrableObjects();
203     nCutoffGroups_ -= mol->getNCutoffGroups();
204     nConstraints_ -= mol->getNConstraintPairs();
205 gezelter 2
206 gezelter 1287 removeInteractionPairs(mol);
207 gezelter 507 molecules_.erase(mol->getGlobalIndex());
208 gezelter 2
209 gezelter 507 delete mol;
210 gezelter 246
211 gezelter 507 return true;
212 gezelter 246 } else {
213 gezelter 507 return false;
214 gezelter 246 }
215 gezelter 507 }
216 gezelter 246
217    
218 gezelter 507 Molecule* SimInfo::beginMolecule(MoleculeIterator& i) {
219 gezelter 246 i = molecules_.begin();
220     return i == molecules_.end() ? NULL : i->second;
221 gezelter 507 }
222 gezelter 246
223 gezelter 507 Molecule* SimInfo::nextMolecule(MoleculeIterator& i) {
224 gezelter 246 ++i;
225     return i == molecules_.end() ? NULL : i->second;
226 gezelter 507 }
227 gezelter 2
228    
229 gezelter 507 void SimInfo::calcNdf() {
230 gezelter 246 int ndf_local;
231     MoleculeIterator i;
232 gezelter 1528 vector<StuntDouble*>::iterator j;
233 gezelter 246 Molecule* mol;
234     StuntDouble* integrableObject;
235 gezelter 2
236 gezelter 246 ndf_local = 0;
237    
238     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
239 gezelter 507 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
240     integrableObject = mol->nextIntegrableObject(j)) {
241 gezelter 2
242 gezelter 507 ndf_local += 3;
243 gezelter 2
244 gezelter 507 if (integrableObject->isDirectional()) {
245     if (integrableObject->isLinear()) {
246     ndf_local += 2;
247     } else {
248     ndf_local += 3;
249     }
250     }
251 gezelter 246
252 tim 770 }
253     }
254 gezelter 246
255     // n_constraints is local, so subtract them on each processor
256     ndf_local -= nConstraints_;
257    
258     #ifdef IS_MPI
259     MPI_Allreduce(&ndf_local,&ndf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
260     #else
261     ndf_ = ndf_local;
262     #endif
263    
264     // nZconstraints_ is global, as are the 3 COM translations for the
265     // entire system:
266     ndf_ = ndf_ - 3 - nZconstraint_;
267    
268 gezelter 507 }
269 gezelter 2
270 gezelter 945 int SimInfo::getFdf() {
271     #ifdef IS_MPI
272     MPI_Allreduce(&fdf_local,&fdf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
273     #else
274     fdf_ = fdf_local;
275     #endif
276     return fdf_;
277     }
278    
279 gezelter 507 void SimInfo::calcNdfRaw() {
280 gezelter 246 int ndfRaw_local;
281 gezelter 2
282 gezelter 246 MoleculeIterator i;
283 gezelter 1528 vector<StuntDouble*>::iterator j;
284 gezelter 246 Molecule* mol;
285     StuntDouble* integrableObject;
286    
287     // Raw degrees of freedom that we have to set
288     ndfRaw_local = 0;
289    
290     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
291 gezelter 507 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
292     integrableObject = mol->nextIntegrableObject(j)) {
293 gezelter 246
294 gezelter 507 ndfRaw_local += 3;
295 gezelter 246
296 gezelter 507 if (integrableObject->isDirectional()) {
297     if (integrableObject->isLinear()) {
298     ndfRaw_local += 2;
299     } else {
300     ndfRaw_local += 3;
301     }
302     }
303 gezelter 246
304 gezelter 507 }
305 gezelter 246 }
306    
307     #ifdef IS_MPI
308     MPI_Allreduce(&ndfRaw_local,&ndfRaw_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
309     #else
310     ndfRaw_ = ndfRaw_local;
311     #endif
312 gezelter 507 }
313 gezelter 2
314 gezelter 507 void SimInfo::calcNdfTrans() {
315 gezelter 246 int ndfTrans_local;
316 gezelter 2
317 gezelter 246 ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_;
318 gezelter 2
319    
320 gezelter 246 #ifdef IS_MPI
321     MPI_Allreduce(&ndfTrans_local,&ndfTrans_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
322     #else
323     ndfTrans_ = ndfTrans_local;
324     #endif
325 gezelter 2
326 gezelter 246 ndfTrans_ = ndfTrans_ - 3 - nZconstraint_;
327    
328 gezelter 507 }
329 gezelter 2
330 gezelter 1287 void SimInfo::addInteractionPairs(Molecule* mol) {
331     ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
332 gezelter 1528 vector<Bond*>::iterator bondIter;
333     vector<Bend*>::iterator bendIter;
334     vector<Torsion*>::iterator torsionIter;
335     vector<Inversion*>::iterator inversionIter;
336 gezelter 246 Bond* bond;
337     Bend* bend;
338     Torsion* torsion;
339 gezelter 1277 Inversion* inversion;
340 gezelter 246 int a;
341     int b;
342     int c;
343     int d;
344 tim 749
345 gezelter 1287 // atomGroups can be used to add special interaction maps between
346     // groups of atoms that are in two separate rigid bodies.
347     // However, most site-site interactions between two rigid bodies
348     // are probably not special, just the ones between the physically
349     // bonded atoms. Interactions *within* a single rigid body should
350     // always be excluded. These are done at the bottom of this
351     // function.
352    
353 gezelter 1528 map<int, set<int> > atomGroups;
354 tim 749 Molecule::RigidBodyIterator rbIter;
355     RigidBody* rb;
356     Molecule::IntegrableObjectIterator ii;
357     StuntDouble* integrableObject;
358 gezelter 246
359 gezelter 1287 for (integrableObject = mol->beginIntegrableObject(ii);
360     integrableObject != NULL;
361     integrableObject = mol->nextIntegrableObject(ii)) {
362    
363 tim 749 if (integrableObject->isRigidBody()) {
364 gezelter 1287 rb = static_cast<RigidBody*>(integrableObject);
365 gezelter 1528 vector<Atom*> atoms = rb->getAtoms();
366     set<int> rigidAtoms;
367 gezelter 1287 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
368     rigidAtoms.insert(atoms[i]->getGlobalIndex());
369     }
370     for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
371 gezelter 1528 atomGroups.insert(map<int, set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
372 gezelter 1287 }
373 tim 749 } else {
374 gezelter 1528 set<int> oneAtomSet;
375 tim 749 oneAtomSet.insert(integrableObject->getGlobalIndex());
376 gezelter 1528 atomGroups.insert(map<int, set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));
377 tim 749 }
378     }
379 gezelter 1287
380     for (bond= mol->beginBond(bondIter); bond != NULL;
381     bond = mol->nextBond(bondIter)) {
382 tim 749
383 gezelter 1287 a = bond->getAtomA()->getGlobalIndex();
384     b = bond->getAtomB()->getGlobalIndex();
385 tim 749
386 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
387     oneTwoInteractions_.addPair(a, b);
388     } else {
389     excludedInteractions_.addPair(a, b);
390     }
391 gezelter 246 }
392 gezelter 2
393 gezelter 1287 for (bend= mol->beginBend(bendIter); bend != NULL;
394     bend = mol->nextBend(bendIter)) {
395    
396 gezelter 507 a = bend->getAtomA()->getGlobalIndex();
397     b = bend->getAtomB()->getGlobalIndex();
398     c = bend->getAtomC()->getGlobalIndex();
399 gezelter 1287
400     if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
401     oneTwoInteractions_.addPair(a, b);
402     oneTwoInteractions_.addPair(b, c);
403     } else {
404     excludedInteractions_.addPair(a, b);
405     excludedInteractions_.addPair(b, c);
406     }
407 gezelter 2
408 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
409     oneThreeInteractions_.addPair(a, c);
410     } else {
411     excludedInteractions_.addPair(a, c);
412     }
413 gezelter 246 }
414 gezelter 2
415 gezelter 1287 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL;
416     torsion = mol->nextTorsion(torsionIter)) {
417    
418 gezelter 507 a = torsion->getAtomA()->getGlobalIndex();
419     b = torsion->getAtomB()->getGlobalIndex();
420     c = torsion->getAtomC()->getGlobalIndex();
421 gezelter 1287 d = torsion->getAtomD()->getGlobalIndex();
422 cli2 1290
423 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
424     oneTwoInteractions_.addPair(a, b);
425     oneTwoInteractions_.addPair(b, c);
426     oneTwoInteractions_.addPair(c, d);
427     } else {
428     excludedInteractions_.addPair(a, b);
429     excludedInteractions_.addPair(b, c);
430     excludedInteractions_.addPair(c, d);
431     }
432 gezelter 2
433 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
434     oneThreeInteractions_.addPair(a, c);
435     oneThreeInteractions_.addPair(b, d);
436     } else {
437     excludedInteractions_.addPair(a, c);
438     excludedInteractions_.addPair(b, d);
439     }
440 tim 749
441 gezelter 1287 if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
442     oneFourInteractions_.addPair(a, d);
443     } else {
444     excludedInteractions_.addPair(a, d);
445     }
446 gezelter 2 }
447    
448 gezelter 1277 for (inversion= mol->beginInversion(inversionIter); inversion != NULL;
449     inversion = mol->nextInversion(inversionIter)) {
450 gezelter 1287
451 gezelter 1277 a = inversion->getAtomA()->getGlobalIndex();
452     b = inversion->getAtomB()->getGlobalIndex();
453     c = inversion->getAtomC()->getGlobalIndex();
454     d = inversion->getAtomD()->getGlobalIndex();
455    
456 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
457     oneTwoInteractions_.addPair(a, b);
458     oneTwoInteractions_.addPair(a, c);
459     oneTwoInteractions_.addPair(a, d);
460     } else {
461     excludedInteractions_.addPair(a, b);
462     excludedInteractions_.addPair(a, c);
463     excludedInteractions_.addPair(a, d);
464     }
465 gezelter 1277
466 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
467     oneThreeInteractions_.addPair(b, c);
468     oneThreeInteractions_.addPair(b, d);
469     oneThreeInteractions_.addPair(c, d);
470     } else {
471     excludedInteractions_.addPair(b, c);
472     excludedInteractions_.addPair(b, d);
473     excludedInteractions_.addPair(c, d);
474     }
475 gezelter 1277 }
476    
477 gezelter 1287 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
478     rb = mol->nextRigidBody(rbIter)) {
479 gezelter 1528 vector<Atom*> atoms = rb->getAtoms();
480 gezelter 1287 for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) {
481     for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
482 gezelter 507 a = atoms[i]->getGlobalIndex();
483     b = atoms[j]->getGlobalIndex();
484 gezelter 1287 excludedInteractions_.addPair(a, b);
485 gezelter 507 }
486     }
487 tim 430 }
488    
489 gezelter 507 }
490 gezelter 246
491 gezelter 1287 void SimInfo::removeInteractionPairs(Molecule* mol) {
492     ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
493 gezelter 1528 vector<Bond*>::iterator bondIter;
494     vector<Bend*>::iterator bendIter;
495     vector<Torsion*>::iterator torsionIter;
496     vector<Inversion*>::iterator inversionIter;
497 gezelter 246 Bond* bond;
498     Bend* bend;
499     Torsion* torsion;
500 gezelter 1277 Inversion* inversion;
501 gezelter 246 int a;
502     int b;
503     int c;
504     int d;
505 tim 749
506 gezelter 1528 map<int, set<int> > atomGroups;
507 tim 749 Molecule::RigidBodyIterator rbIter;
508     RigidBody* rb;
509     Molecule::IntegrableObjectIterator ii;
510     StuntDouble* integrableObject;
511 gezelter 246
512 gezelter 1287 for (integrableObject = mol->beginIntegrableObject(ii);
513     integrableObject != NULL;
514     integrableObject = mol->nextIntegrableObject(ii)) {
515    
516 tim 749 if (integrableObject->isRigidBody()) {
517 gezelter 1287 rb = static_cast<RigidBody*>(integrableObject);
518 gezelter 1528 vector<Atom*> atoms = rb->getAtoms();
519     set<int> rigidAtoms;
520 gezelter 1287 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
521     rigidAtoms.insert(atoms[i]->getGlobalIndex());
522     }
523     for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
524 gezelter 1528 atomGroups.insert(map<int, set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
525 gezelter 1287 }
526 tim 749 } else {
527 gezelter 1528 set<int> oneAtomSet;
528 tim 749 oneAtomSet.insert(integrableObject->getGlobalIndex());
529 gezelter 1528 atomGroups.insert(map<int, set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));
530 tim 749 }
531     }
532    
533 gezelter 1287 for (bond= mol->beginBond(bondIter); bond != NULL;
534     bond = mol->nextBond(bondIter)) {
535    
536     a = bond->getAtomA()->getGlobalIndex();
537     b = bond->getAtomB()->getGlobalIndex();
538 tim 749
539 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
540     oneTwoInteractions_.removePair(a, b);
541     } else {
542     excludedInteractions_.removePair(a, b);
543     }
544 gezelter 2 }
545 gezelter 246
546 gezelter 1287 for (bend= mol->beginBend(bendIter); bend != NULL;
547     bend = mol->nextBend(bendIter)) {
548    
549 gezelter 507 a = bend->getAtomA()->getGlobalIndex();
550     b = bend->getAtomB()->getGlobalIndex();
551     c = bend->getAtomC()->getGlobalIndex();
552 gezelter 1287
553     if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
554     oneTwoInteractions_.removePair(a, b);
555     oneTwoInteractions_.removePair(b, c);
556     } else {
557     excludedInteractions_.removePair(a, b);
558     excludedInteractions_.removePair(b, c);
559     }
560 gezelter 246
561 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
562     oneThreeInteractions_.removePair(a, c);
563     } else {
564     excludedInteractions_.removePair(a, c);
565     }
566 gezelter 2 }
567 gezelter 246
568 gezelter 1287 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL;
569     torsion = mol->nextTorsion(torsionIter)) {
570    
571 gezelter 507 a = torsion->getAtomA()->getGlobalIndex();
572     b = torsion->getAtomB()->getGlobalIndex();
573     c = torsion->getAtomC()->getGlobalIndex();
574 gezelter 1287 d = torsion->getAtomD()->getGlobalIndex();
575    
576     if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
577     oneTwoInteractions_.removePair(a, b);
578     oneTwoInteractions_.removePair(b, c);
579     oneTwoInteractions_.removePair(c, d);
580     } else {
581     excludedInteractions_.removePair(a, b);
582     excludedInteractions_.removePair(b, c);
583     excludedInteractions_.removePair(c, d);
584     }
585 gezelter 246
586 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
587     oneThreeInteractions_.removePair(a, c);
588     oneThreeInteractions_.removePair(b, d);
589     } else {
590     excludedInteractions_.removePair(a, c);
591     excludedInteractions_.removePair(b, d);
592     }
593 tim 749
594 gezelter 1287 if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
595     oneFourInteractions_.removePair(a, d);
596     } else {
597     excludedInteractions_.removePair(a, d);
598     }
599     }
600 tim 749
601 gezelter 1287 for (inversion= mol->beginInversion(inversionIter); inversion != NULL;
602     inversion = mol->nextInversion(inversionIter)) {
603 tim 749
604 gezelter 1277 a = inversion->getAtomA()->getGlobalIndex();
605     b = inversion->getAtomB()->getGlobalIndex();
606     c = inversion->getAtomC()->getGlobalIndex();
607     d = inversion->getAtomD()->getGlobalIndex();
608    
609 gezelter 1287 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
610     oneTwoInteractions_.removePair(a, b);
611     oneTwoInteractions_.removePair(a, c);
612     oneTwoInteractions_.removePair(a, d);
613     } else {
614     excludedInteractions_.removePair(a, b);
615     excludedInteractions_.removePair(a, c);
616     excludedInteractions_.removePair(a, d);
617     }
618 gezelter 1277
619 gezelter 1287 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
620     oneThreeInteractions_.removePair(b, c);
621     oneThreeInteractions_.removePair(b, d);
622     oneThreeInteractions_.removePair(c, d);
623     } else {
624     excludedInteractions_.removePair(b, c);
625     excludedInteractions_.removePair(b, d);
626     excludedInteractions_.removePair(c, d);
627     }
628 gezelter 1277 }
629    
630 gezelter 1287 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
631     rb = mol->nextRigidBody(rbIter)) {
632 gezelter 1528 vector<Atom*> atoms = rb->getAtoms();
633 gezelter 1287 for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) {
634     for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
635 gezelter 507 a = atoms[i]->getGlobalIndex();
636     b = atoms[j]->getGlobalIndex();
637 gezelter 1287 excludedInteractions_.removePair(a, b);
638 gezelter 507 }
639     }
640 tim 430 }
641 gezelter 1287
642 gezelter 507 }
643 gezelter 1287
644    
645 gezelter 507 void SimInfo::addMoleculeStamp(MoleculeStamp* molStamp, int nmol) {
646 gezelter 246 int curStampId;
647 gezelter 1287
648 gezelter 246 //index from 0
649     curStampId = moleculeStamps_.size();
650 gezelter 2
651 gezelter 246 moleculeStamps_.push_back(molStamp);
652     molStampIds_.insert(molStampIds_.end(), nmol, curStampId);
653 gezelter 507 }
654 gezelter 2
655 gezelter 1530
656     /**
657     * update
658     *
659     * Performs the global checks and variable settings after the objects have been
660     * created.
661     *
662     */
663 gezelter 507 void SimInfo::update() {
664 gezelter 1530
665     setupSimVariables();
666     setupCutoffs();
667     setupSwitching();
668     setupElectrostatics();
669     setupNeighborlists();
670 gezelter 2
671 gezelter 246 #ifdef IS_MPI
672     setupFortranParallel();
673     #endif
674     setupFortranSim();
675 gezelter 1528 fortranInitialized_ = true;
676 gezelter 2
677 gezelter 246 calcNdf();
678     calcNdfRaw();
679     calcNdfTrans();
680 gezelter 507 }
681 gezelter 1528
682     set<AtomType*> SimInfo::getSimulatedAtomTypes() {
683 gezelter 246 SimInfo::MoleculeIterator mi;
684     Molecule* mol;
685     Molecule::AtomIterator ai;
686     Atom* atom;
687 gezelter 1528 set<AtomType*> atomTypes;
688    
689 gezelter 1529 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
690 gezelter 507 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
691     atomTypes.insert(atom->getAtomType());
692 gezelter 1529 }
693     }
694 gezelter 246 return atomTypes;
695 gezelter 507 }
696 gezelter 2
697 gezelter 1528 /**
698 gezelter 1530 * setupCutoffs
699 gezelter 1528 *
700 gezelter 1530 * Sets the values of cutoffRadius and cutoffMethod
701     *
702     * cutoffRadius : realType
703 gezelter 1528 * If the cutoffRadius was explicitly set, use that value.
704     * If the cutoffRadius was not explicitly set:
705     * Are there electrostatic atoms? Use 12.0 Angstroms.
706     * No electrostatic atoms? Poll the atom types present in the
707     * simulation for suggested cutoff values (e.g. 2.5 * sigma).
708     * Use the maximum suggested value that was found.
709 gezelter 1530 *
710     * cutoffMethod : (one of HARD, SWITCHED, SHIFTED_FORCE, SHIFTED_POTENTIAL)
711     * If cutoffMethod was explicitly set, use that choice.
712     * If cutoffMethod was not explicitly set, use SHIFTED_FORCE
713 gezelter 1528 */
714 gezelter 1530 void SimInfo::setupCutoffs() {
715 gezelter 2
716 gezelter 1528 if (simParams_->haveCutoffRadius()) {
717     cutoffRadius_ = simParams_->getCutoffRadius();
718     } else {
719     if (usesElectrostaticAtoms_) {
720     sprintf(painCave.errMsg,
721 gezelter 1530 "SimInfo: No value was set for the cutoffRadius.\n"
722 gezelter 1528 "\tOpenMD will use a default value of 12.0 angstroms"
723     "\tfor the cutoffRadius.\n");
724     painCave.isFatal = 0;
725 gezelter 1530 painCave.severity = OPENMD_INFO;
726 gezelter 1528 simError();
727     cutoffRadius_ = 12.0;
728     } else {
729     RealType thisCut;
730     set<AtomType*>::iterator i;
731     set<AtomType*> atomTypes;
732     atomTypes = getSimulatedAtomTypes();
733     for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
734     thisCut = InteractionManager::Instance()->getSuggestedCutoffRadius((*i));
735     cutoffRadius_ = max(thisCut, cutoffRadius_);
736     }
737     sprintf(painCave.errMsg,
738 gezelter 1530 "SimInfo: No value was set for the cutoffRadius.\n"
739 gezelter 1528 "\tOpenMD will use %lf angstroms.\n",
740     cutoffRadius_);
741     painCave.isFatal = 0;
742 gezelter 1530 painCave.severity = OPENMD_INFO;
743 gezelter 1528 simError();
744     }
745     }
746 gezelter 1126
747 gezelter 1530 map<string, CutoffMethod> stringToCutoffMethod;
748     stringToCutoffMethod["HARD"] = HARD;
749     stringToCutoffMethod["SWITCHING_FUNCTION"] = SWITCHING_FUNCTION;
750     stringToCutoffMethod["SHIFTED_POTENTIAL"] = SHIFTED_POTENTIAL;
751     stringToCutoffMethod["SHIFTED_FORCE"] = SHIFTED_FORCE;
752    
753     if (simParams_->haveCutoffMethod()) {
754     string cutMeth = toUpperCopy(simParams_->getCutoffMethod());
755     map<string, CutoffMethod>::iterator i;
756     i = stringToCutoffMethod.find(cutMeth);
757     if (i == stringToCutoffMethod.end()) {
758     sprintf(painCave.errMsg,
759     "SimInfo: Could not find chosen cutoffMethod %s\n"
760     "\tShould be one of: "
761     "HARD, SWITCHING_FUNCTION, SHIFTED_POTENTIAL, or SHIFTED_FORCE\n",
762     cutMeth.c_str());
763     painCave.isFatal = 1;
764     painCave.severity = OPENMD_ERROR;
765     simError();
766     } else {
767     cutoffMethod_ = i->second;
768     }
769     } else {
770     sprintf(painCave.errMsg,
771     "SimInfo: No value was set for the cutoffMethod.\n"
772     "\tOpenMD will use SHIFTED_FORCE.\n");
773     painCave.isFatal = 0;
774     painCave.severity = OPENMD_INFO;
775     simError();
776     cutoffMethod_ = SHIFTED_FORCE;
777     }
778 gezelter 1528 }
779    
780     /**
781 gezelter 1530 * setupSwitching
782 gezelter 1528 *
783 gezelter 1530 * Sets the values of switchingRadius and
784 gezelter 1528 * If the switchingRadius was explicitly set, use that value (but check it)
785     * If the switchingRadius was not explicitly set: use 0.85 * cutoffRadius_
786     */
787 gezelter 1530 void SimInfo::setupSwitching() {
788 gezelter 1528
789     if (simParams_->haveSwitchingRadius()) {
790     switchingRadius_ = simParams_->getSwitchingRadius();
791     if (switchingRadius_ > cutoffRadius_) {
792     sprintf(painCave.errMsg,
793 gezelter 1530 "SimInfo: switchingRadius (%f) is larger than cutoffRadius(%f)\n",
794 gezelter 1528 switchingRadius_, cutoffRadius_);
795     painCave.isFatal = 1;
796 gezelter 1530 painCave.severity = OPENMD_ERROR;
797 gezelter 1528 simError();
798 chrisfen 691 }
799 gezelter 1528 } else {
800     switchingRadius_ = 0.85 * cutoffRadius_;
801     sprintf(painCave.errMsg,
802 gezelter 1530 "SimInfo: No value was set for the switchingRadius.\n"
803 gezelter 1528 "\tOpenMD will use a default value of 85 percent of the cutoffRadius.\n"
804     "\tswitchingRadius = %f. for this simulation\n", switchingRadius_);
805     painCave.isFatal = 0;
806 gezelter 1530 painCave.severity = OPENMD_WARNING;
807 gezelter 1528 simError();
808 gezelter 1530 }
809    
810     if (simParams_->haveSwitchingFunctionType()) {
811     string funcType = simParams_->getSwitchingFunctionType();
812     toUpper(funcType);
813     if (funcType == "CUBIC") {
814 gezelter 1534 sft_ = cubic;
815 gezelter 1530 } else {
816     if (funcType == "FIFTH_ORDER_POLYNOMIAL") {
817 gezelter 1534 sft_ = fifth_order_poly;
818 gezelter 1530 } else {
819     // throw error
820     sprintf( painCave.errMsg,
821     "SimInfo : Unknown switchingFunctionType. (Input file specified %s .)\n"
822     "\tswitchingFunctionType must be one of: "
823     "\"cubic\" or \"fifth_order_polynomial\".",
824     funcType.c_str() );
825     painCave.isFatal = 1;
826     painCave.severity = OPENMD_ERROR;
827     simError();
828     }
829     }
830     }
831 gezelter 1528 }
832 chrisfen 611
833 gezelter 1528 /**
834 gezelter 1534 * setupNeighborlists
835 gezelter 1528 *
836     * If the skinThickness was explicitly set, use that value (but check it)
837     * If the skinThickness was not explicitly set: use 1.0 angstroms
838     */
839 gezelter 1534 void SimInfo::setupNeighborlists() {
840 gezelter 1528 if (simParams_->haveSkinThickness()) {
841     skinThickness_ = simParams_->getSkinThickness();
842     } else {
843     skinThickness_ = 1.0;
844     sprintf(painCave.errMsg,
845 gezelter 1534 "SimInfo: No value was set for the skinThickness.\n"
846 gezelter 1528 "\tOpenMD will use a default value of %f Angstroms\n"
847     "\tfor this simulation\n", skinThickness_);
848 gezelter 1534 painCave.severity = OPENMD_INFO;
849 gezelter 1528 painCave.isFatal = 0;
850     simError();
851     }
852     }
853    
854 gezelter 1534 void SimInfo::setupSimVariables() {
855     useAtomicVirial_ = simParams_->getUseAtomicVirial();
856     // we only call setAccumulateBoxDipole if the accumulateBoxDipole parameter is true
857     calcBoxDipole_ = false;
858     if ( simParams_->haveAccumulateBoxDipole() )
859     if ( simParams_->getAccumulateBoxDipole() ) {
860     calcBoxDipole_ = true;
861     }
862    
863 gezelter 1528 set<AtomType*>::iterator i;
864     set<AtomType*> atomTypes;
865 gezelter 1534 atomTypes = getSimulatedAtomTypes();
866 gezelter 1528 int usesElectrostatic = 0;
867     int usesMetallic = 0;
868     int usesDirectional = 0;
869 gezelter 246 //loop over all of the atom types
870     for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
871 gezelter 1528 usesElectrostatic |= (*i)->isElectrostatic();
872     usesMetallic |= (*i)->isMetal();
873     usesDirectional |= (*i)->isDirectional();
874 gezelter 246 }
875 gezelter 2
876 gezelter 246 #ifdef IS_MPI
877     int temp;
878 gezelter 1528 temp = usesDirectional;
879     MPI_Allreduce(&temp, &usesDirectionalAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
880 gezelter 2
881 gezelter 1528 temp = usesMetallic;
882     MPI_Allreduce(&temp, &usesMetallicAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
883 gezelter 2
884 gezelter 1528 temp = usesElectrostatic;
885     MPI_Allreduce(&temp, &usesElectrostaticAtoms_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
886 gezelter 2 #endif
887 gezelter 1528 fInfo_.SIM_uses_PBC = usesPeriodicBoundaries_;
888     fInfo_.SIM_uses_DirectionalAtoms = usesDirectionalAtoms_;
889     fInfo_.SIM_uses_MetallicAtoms = usesMetallicAtoms_;
890     fInfo_.SIM_requires_SkipCorrection = usesElectrostaticAtoms_;
891     fInfo_.SIM_requires_SelfCorrection = usesElectrostaticAtoms_;
892     fInfo_.SIM_uses_AtomicVirial = usesAtomicVirial_;
893 gezelter 507 }
894 gezelter 2
895 gezelter 507 void SimInfo::setupFortranSim() {
896 gezelter 246 int isError;
897 gezelter 1287 int nExclude, nOneTwo, nOneThree, nOneFour;
898 gezelter 1528 vector<int> fortranGlobalGroupMembership;
899 gezelter 246
900 gezelter 1528 notifyFortranSkinThickness(&skinThickness_);
901    
902     int ljsp = cutoffMethod_ == SHIFTED_POTENTIAL ? 1 : 0;
903     int ljsf = cutoffMethod_ == SHIFTED_FORCE ? 1 : 0;
904     notifyFortranCutoffs(&cutoffRadius_, &switchingRadius_, &ljsp, &ljsf);
905    
906 gezelter 246 isError = 0;
907 gezelter 2
908 gezelter 246 //globalGroupMembership_ is filled by SimCreator
909     for (int i = 0; i < nGlobalAtoms_; i++) {
910 gezelter 507 fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1);
911 gezelter 246 }
912 gezelter 2
913 gezelter 246 //calculate mass ratio of cutoff group
914 gezelter 1528 vector<RealType> mfact;
915 gezelter 246 SimInfo::MoleculeIterator mi;
916     Molecule* mol;
917     Molecule::CutoffGroupIterator ci;
918     CutoffGroup* cg;
919     Molecule::AtomIterator ai;
920     Atom* atom;
921 tim 963 RealType totalMass;
922 gezelter 246
923     //to avoid memory reallocation, reserve enough space for mfact
924     mfact.reserve(getNCutoffGroups());
925 gezelter 2
926 gezelter 246 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
927 gezelter 507 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
928 gezelter 2
929 gezelter 507 totalMass = cg->getMass();
930     for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
931 chrisfen 645 // Check for massless groups - set mfact to 1 if true
932     if (totalMass != 0)
933     mfact.push_back(atom->getMass()/totalMass);
934     else
935     mfact.push_back( 1.0 );
936 gezelter 507 }
937     }
938 gezelter 246 }
939 gezelter 2
940 gezelter 246 //fill ident array of local atoms (it is actually ident of AtomType, it is so confusing !!!)
941 gezelter 1528 vector<int> identArray;
942 gezelter 2
943 gezelter 246 //to avoid memory reallocation, reserve enough space identArray
944     identArray.reserve(getNAtoms());
945    
946     for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
947 gezelter 507 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
948     identArray.push_back(atom->getIdent());
949     }
950 gezelter 246 }
951 gezelter 2
952 gezelter 246 //fill molMembershipArray
953     //molMembershipArray is filled by SimCreator
954 gezelter 1528 vector<int> molMembershipArray(nGlobalAtoms_);
955 gezelter 246 for (int i = 0; i < nGlobalAtoms_; i++) {
956 gezelter 507 molMembershipArray[i] = globalMolMembership_[i] + 1;
957 gezelter 246 }
958    
959     //setup fortran simulation
960 gezelter 1287
961     nExclude = excludedInteractions_.getSize();
962     nOneTwo = oneTwoInteractions_.getSize();
963     nOneThree = oneThreeInteractions_.getSize();
964     nOneFour = oneFourInteractions_.getSize();
965    
966     int* excludeList = excludedInteractions_.getPairList();
967     int* oneTwoList = oneTwoInteractions_.getPairList();
968     int* oneThreeList = oneThreeInteractions_.getPairList();
969     int* oneFourList = oneFourInteractions_.getPairList();
970    
971 gezelter 1241 setFortranSim( &fInfo_, &nGlobalAtoms_, &nAtoms_, &identArray[0],
972 gezelter 1287 &nExclude, excludeList,
973     &nOneTwo, oneTwoList,
974     &nOneThree, oneThreeList,
975     &nOneFour, oneFourList,
976 gezelter 1241 &molMembershipArray[0], &mfact[0], &nCutoffGroups_,
977     &fortranGlobalGroupMembership[0], &isError);
978    
979 gezelter 246 if( isError ){
980 gezelter 1241
981 gezelter 507 sprintf( painCave.errMsg,
982     "There was an error setting the simulation information in fortran.\n" );
983     painCave.isFatal = 1;
984 gezelter 1390 painCave.severity = OPENMD_ERROR;
985 gezelter 507 simError();
986 gezelter 246 }
987 gezelter 1241
988    
989 gezelter 246 sprintf( checkPointMsg,
990 gezelter 507 "succesfully sent the simulation information to fortran.\n");
991 gezelter 1241
992     errorCheckPoint();
993    
994 chuckv 1095 // Setup number of neighbors in neighbor list if present
995     if (simParams_->haveNeighborListNeighbors()) {
996 chuckv 1121 int nlistNeighbors = simParams_->getNeighborListNeighbors();
997     setNeighbors(&nlistNeighbors);
998 chuckv 1095 }
999    
1000    
1001 gezelter 507 }
1002 gezelter 2
1003    
1004 gezelter 507 void SimInfo::setupFortranParallel() {
1005 gezelter 1241 #ifdef IS_MPI
1006 gezelter 246 //SimInfo is responsible for creating localToGlobalAtomIndex and localToGlobalGroupIndex
1007 gezelter 1528 vector<int> localToGlobalAtomIndex(getNAtoms(), 0);
1008     vector<int> localToGlobalCutoffGroupIndex;
1009 gezelter 246 SimInfo::MoleculeIterator mi;
1010     Molecule::AtomIterator ai;
1011     Molecule::CutoffGroupIterator ci;
1012     Molecule* mol;
1013     Atom* atom;
1014     CutoffGroup* cg;
1015     mpiSimData parallelData;
1016     int isError;
1017 gezelter 2
1018 gezelter 246 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
1019 gezelter 2
1020 gezelter 507 //local index(index in DataStorge) of atom is important
1021     for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
1022     localToGlobalAtomIndex[atom->getLocalIndex()] = atom->getGlobalIndex() + 1;
1023     }
1024 gezelter 2
1025 gezelter 507 //local index of cutoff group is trivial, it only depends on the order of travesing
1026     for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
1027     localToGlobalCutoffGroupIndex.push_back(cg->getGlobalIndex() + 1);
1028     }
1029 gezelter 246
1030     }
1031 gezelter 2
1032 gezelter 246 //fill up mpiSimData struct
1033     parallelData.nMolGlobal = getNGlobalMolecules();
1034     parallelData.nMolLocal = getNMolecules();
1035     parallelData.nAtomsGlobal = getNGlobalAtoms();
1036     parallelData.nAtomsLocal = getNAtoms();
1037     parallelData.nGroupsGlobal = getNGlobalCutoffGroups();
1038     parallelData.nGroupsLocal = getNCutoffGroups();
1039     parallelData.myNode = worldRank;
1040     MPI_Comm_size(MPI_COMM_WORLD, &(parallelData.nProcessors));
1041 gezelter 2
1042 gezelter 246 //pass mpiSimData struct and index arrays to fortran
1043     setFsimParallel(&parallelData, &(parallelData.nAtomsLocal),
1044     &localToGlobalAtomIndex[0], &(parallelData.nGroupsLocal),
1045     &localToGlobalCutoffGroupIndex[0], &isError);
1046 gezelter 2
1047 gezelter 246 if (isError) {
1048 gezelter 507 sprintf(painCave.errMsg,
1049     "mpiRefresh errror: fortran didn't like something we gave it.\n");
1050     painCave.isFatal = 1;
1051     simError();
1052 gezelter 246 }
1053 gezelter 2
1054 gezelter 246 sprintf(checkPointMsg, " mpiRefresh successful.\n");
1055 gezelter 1241 errorCheckPoint();
1056 gezelter 2
1057 gezelter 1241 #endif
1058 gezelter 507 }
1059 chrisfen 143
1060 chuckv 834
1061 chrisfen 998 void SimInfo::setupAccumulateBoxDipole() {
1062    
1063    
1064     }
1065    
1066 gezelter 507 void SimInfo::addProperty(GenericData* genData) {
1067 gezelter 246 properties_.addProperty(genData);
1068 gezelter 507 }
1069 gezelter 2
1070 gezelter 1528 void SimInfo::removeProperty(const string& propName) {
1071 gezelter 246 properties_.removeProperty(propName);
1072 gezelter 507 }
1073 gezelter 2
1074 gezelter 507 void SimInfo::clearProperties() {
1075 gezelter 246 properties_.clearProperties();
1076 gezelter 507 }
1077 gezelter 2
1078 gezelter 1528 vector<string> SimInfo::getPropertyNames() {
1079 gezelter 246 return properties_.getPropertyNames();
1080 gezelter 507 }
1081 gezelter 246
1082 gezelter 1528 vector<GenericData*> SimInfo::getProperties() {
1083 gezelter 246 return properties_.getProperties();
1084 gezelter 507 }
1085 gezelter 2
1086 gezelter 1528 GenericData* SimInfo::getPropertyByName(const string& propName) {
1087 gezelter 246 return properties_.getPropertyByName(propName);
1088 gezelter 507 }
1089 gezelter 2
1090 gezelter 507 void SimInfo::setSnapshotManager(SnapshotManager* sman) {
1091 tim 432 if (sman_ == sman) {
1092 gezelter 507 return;
1093 tim 432 }
1094     delete sman_;
1095 gezelter 246 sman_ = sman;
1096 gezelter 2
1097 gezelter 246 Molecule* mol;
1098     RigidBody* rb;
1099     Atom* atom;
1100     SimInfo::MoleculeIterator mi;
1101     Molecule::RigidBodyIterator rbIter;
1102     Molecule::AtomIterator atomIter;;
1103    
1104     for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
1105    
1106 gezelter 507 for (atom = mol->beginAtom(atomIter); atom != NULL; atom = mol->nextAtom(atomIter)) {
1107     atom->setSnapshotManager(sman_);
1108     }
1109 gezelter 246
1110 gezelter 507 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
1111     rb->setSnapshotManager(sman_);
1112     }
1113 gezelter 246 }
1114 gezelter 2
1115 gezelter 507 }
1116 gezelter 2
1117 gezelter 507 Vector3d SimInfo::getComVel(){
1118 gezelter 246 SimInfo::MoleculeIterator i;
1119     Molecule* mol;
1120 gezelter 2
1121 gezelter 246 Vector3d comVel(0.0);
1122 tim 963 RealType totalMass = 0.0;
1123 gezelter 2
1124 gezelter 246
1125     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1126 tim 963 RealType mass = mol->getMass();
1127 gezelter 507 totalMass += mass;
1128     comVel += mass * mol->getComVel();
1129 gezelter 246 }
1130 gezelter 2
1131 gezelter 246 #ifdef IS_MPI
1132 tim 963 RealType tmpMass = totalMass;
1133 gezelter 246 Vector3d tmpComVel(comVel);
1134 tim 963 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1135     MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1136 gezelter 246 #endif
1137    
1138     comVel /= totalMass;
1139    
1140     return comVel;
1141 gezelter 507 }
1142 gezelter 2
1143 gezelter 507 Vector3d SimInfo::getCom(){
1144 gezelter 246 SimInfo::MoleculeIterator i;
1145     Molecule* mol;
1146 gezelter 2
1147 gezelter 246 Vector3d com(0.0);
1148 tim 963 RealType totalMass = 0.0;
1149 gezelter 246
1150     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1151 tim 963 RealType mass = mol->getMass();
1152 gezelter 507 totalMass += mass;
1153     com += mass * mol->getCom();
1154 gezelter 246 }
1155 gezelter 2
1156     #ifdef IS_MPI
1157 tim 963 RealType tmpMass = totalMass;
1158 gezelter 246 Vector3d tmpCom(com);
1159 tim 963 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1160     MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1161 gezelter 2 #endif
1162    
1163 gezelter 246 com /= totalMass;
1164 gezelter 2
1165 gezelter 246 return com;
1166 gezelter 2
1167 gezelter 507 }
1168 gezelter 246
1169 gezelter 1528 ostream& operator <<(ostream& o, SimInfo& info) {
1170 gezelter 246
1171     return o;
1172 gezelter 507 }
1173 chuckv 555
1174    
1175     /*
1176     Returns center of mass and center of mass velocity in one function call.
1177     */
1178    
1179     void SimInfo::getComAll(Vector3d &com, Vector3d &comVel){
1180     SimInfo::MoleculeIterator i;
1181     Molecule* mol;
1182    
1183    
1184 tim 963 RealType totalMass = 0.0;
1185 chuckv 555
1186 gezelter 246
1187 chuckv 555 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1188 tim 963 RealType mass = mol->getMass();
1189 chuckv 555 totalMass += mass;
1190     com += mass * mol->getCom();
1191     comVel += mass * mol->getComVel();
1192     }
1193    
1194     #ifdef IS_MPI
1195 tim 963 RealType tmpMass = totalMass;
1196 chuckv 555 Vector3d tmpCom(com);
1197     Vector3d tmpComVel(comVel);
1198 tim 963 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1199     MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1200     MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1201 chuckv 555 #endif
1202    
1203     com /= totalMass;
1204     comVel /= totalMass;
1205     }
1206    
1207     /*
1208     Return intertia tensor for entire system and angular momentum Vector.
1209 chuckv 557
1210    
1211     [ Ixx -Ixy -Ixz ]
1212 gezelter 1505 J =| -Iyx Iyy -Iyz |
1213 chuckv 557 [ -Izx -Iyz Izz ]
1214 chuckv 555 */
1215    
1216     void SimInfo::getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum){
1217    
1218    
1219 tim 963 RealType xx = 0.0;
1220     RealType yy = 0.0;
1221     RealType zz = 0.0;
1222     RealType xy = 0.0;
1223     RealType xz = 0.0;
1224     RealType yz = 0.0;
1225 chuckv 555 Vector3d com(0.0);
1226     Vector3d comVel(0.0);
1227    
1228     getComAll(com, comVel);
1229    
1230     SimInfo::MoleculeIterator i;
1231     Molecule* mol;
1232    
1233     Vector3d thisq(0.0);
1234     Vector3d thisv(0.0);
1235    
1236 tim 963 RealType thisMass = 0.0;
1237 chuckv 555
1238    
1239    
1240    
1241     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1242    
1243     thisq = mol->getCom()-com;
1244     thisv = mol->getComVel()-comVel;
1245     thisMass = mol->getMass();
1246     // Compute moment of intertia coefficients.
1247     xx += thisq[0]*thisq[0]*thisMass;
1248     yy += thisq[1]*thisq[1]*thisMass;
1249     zz += thisq[2]*thisq[2]*thisMass;
1250    
1251     // compute products of intertia
1252     xy += thisq[0]*thisq[1]*thisMass;
1253     xz += thisq[0]*thisq[2]*thisMass;
1254     yz += thisq[1]*thisq[2]*thisMass;
1255    
1256     angularMomentum += cross( thisq, thisv ) * thisMass;
1257    
1258     }
1259    
1260    
1261     inertiaTensor(0,0) = yy + zz;
1262     inertiaTensor(0,1) = -xy;
1263     inertiaTensor(0,2) = -xz;
1264     inertiaTensor(1,0) = -xy;
1265 chuckv 557 inertiaTensor(1,1) = xx + zz;
1266 chuckv 555 inertiaTensor(1,2) = -yz;
1267     inertiaTensor(2,0) = -xz;
1268     inertiaTensor(2,1) = -yz;
1269     inertiaTensor(2,2) = xx + yy;
1270    
1271     #ifdef IS_MPI
1272     Mat3x3d tmpI(inertiaTensor);
1273     Vector3d tmpAngMom;
1274 tim 963 MPI_Allreduce(tmpI.getArrayPointer(), inertiaTensor.getArrayPointer(),9,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1275     MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1276 chuckv 555 #endif
1277    
1278     return;
1279     }
1280    
1281     //Returns the angular momentum of the system
1282     Vector3d SimInfo::getAngularMomentum(){
1283    
1284     Vector3d com(0.0);
1285     Vector3d comVel(0.0);
1286     Vector3d angularMomentum(0.0);
1287    
1288     getComAll(com,comVel);
1289    
1290     SimInfo::MoleculeIterator i;
1291     Molecule* mol;
1292    
1293 chuckv 557 Vector3d thisr(0.0);
1294     Vector3d thisp(0.0);
1295 chuckv 555
1296 tim 963 RealType thisMass;
1297 chuckv 555
1298     for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1299 chuckv 557 thisMass = mol->getMass();
1300     thisr = mol->getCom()-com;
1301     thisp = (mol->getComVel()-comVel)*thisMass;
1302 chuckv 555
1303 chuckv 557 angularMomentum += cross( thisr, thisp );
1304    
1305 chuckv 555 }
1306    
1307     #ifdef IS_MPI
1308     Vector3d tmpAngMom;
1309 tim 963 MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1310 chuckv 555 #endif
1311    
1312     return angularMomentum;
1313     }
1314    
1315 tim 1024 StuntDouble* SimInfo::getIOIndexToIntegrableObject(int index) {
1316     return IOIndexToIntegrableObject.at(index);
1317     }
1318    
1319 gezelter 1528 void SimInfo::setIOIndexToIntegrableObject(const vector<StuntDouble*>& v) {
1320 tim 1024 IOIndexToIntegrableObject= v;
1321     }
1322    
1323 chuckv 1103 /* Returns the Volume of the simulation based on a ellipsoid with semi-axes
1324     based on the radius of gyration V=4/3*Pi*R_1*R_2*R_3
1325     where R_i are related to the principle inertia moments R_i = sqrt(C*I_i/N), this reduces to
1326     V = 4/3*Pi*(C/N)^3/2*sqrt(det(I)). See S.E. Baltazar et. al. Comp. Mat. Sci. 37 (2006) 526-536.
1327     */
1328     void SimInfo::getGyrationalVolume(RealType &volume){
1329     Mat3x3d intTensor;
1330     RealType det;
1331     Vector3d dummyAngMom;
1332     RealType sysconstants;
1333     RealType geomCnst;
1334    
1335     geomCnst = 3.0/2.0;
1336     /* Get the inertial tensor and angular momentum for free*/
1337     getInertiaTensor(intTensor,dummyAngMom);
1338    
1339     det = intTensor.determinant();
1340     sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_;
1341     volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(det);
1342     return;
1343     }
1344    
1345     void SimInfo::getGyrationalVolume(RealType &volume, RealType &detI){
1346     Mat3x3d intTensor;
1347     Vector3d dummyAngMom;
1348     RealType sysconstants;
1349     RealType geomCnst;
1350    
1351     geomCnst = 3.0/2.0;
1352     /* Get the inertial tensor and angular momentum for free*/
1353     getInertiaTensor(intTensor,dummyAngMom);
1354    
1355     detI = intTensor.determinant();
1356     sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_;
1357     volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(detI);
1358     return;
1359     }
1360 tim 1024 /*
1361 gezelter 1528 void SimInfo::setStuntDoubleFromGlobalIndex(vector<StuntDouble*> v) {
1362 tim 1024 assert( v.size() == nAtoms_ + nRigidBodies_);
1363     sdByGlobalIndex_ = v;
1364     }
1365    
1366     StuntDouble* SimInfo::getStuntDoubleFromGlobalIndex(int index) {
1367     //assert(index < nAtoms_ + nRigidBodies_);
1368     return sdByGlobalIndex_.at(index);
1369     }
1370     */
1371 gezelter 1528 int SimInfo::getNGlobalConstraints() {
1372     int nGlobalConstraints;
1373     #ifdef IS_MPI
1374     MPI_Allreduce(&nConstraints_, &nGlobalConstraints, 1, MPI_INT, MPI_SUM,
1375     MPI_COMM_WORLD);
1376     #else
1377     nGlobalConstraints = nConstraints_;
1378     #endif
1379     return nGlobalConstraints;
1380     }
1381    
1382 gezelter 1390 }//end namespace OpenMD
1383 gezelter 246

Properties

Name Value
svn:keywords Author Id Revision Date