ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimInfo.cpp
Revision: 1550
Committed: Wed Apr 27 21:49:59 2011 UTC (14 years ago) by gezelter
File size: 40694 byte(s)
Log Message:
more fortran removal

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

Properties

Name Value
svn:keywords Author Id Revision Date