ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.cpp
Revision: 1528
Committed: Fri Dec 17 20:11:05 2010 UTC (14 years, 4 months ago) by gezelter
File size: 43082 byte(s)
Log Message:
Doesn't build yet, but getting much closer...


File Contents

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

Properties

Name Value
svn:keywords Author Id Revision Date