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

Properties

Name Value
svn:keywords Author Id Revision Date