ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.cpp
Revision: 1530
Committed: Tue Dec 28 21:47:55 2010 UTC (14 years, 4 months ago) by gezelter
File size: 44905 byte(s)
Log Message:
Moved switching functions and force options over to the C++ side, and
removed them from Fortran.

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

Properties

Name Value
svn:keywords Author Id Revision Date