ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/SimInfo.cpp
Revision: 1287
Committed: Wed Sep 10 18:11:32 2008 UTC (16 years, 7 months ago) by gezelter
File size: 52334 byte(s)
Log Message:
more changes for 1-2, 1-3, 1-4 interactions plus some initialization-ordering
fixes to make gcc -Wall happier.

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. Acknowledgement of the program authors must be made in any
10 * publication of scientific results based in part on use of the
11 * program. An acceptable form of acknowledgement is citation of
12 * the article in which the program was described (Matthew
13 * A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14 * J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15 * Parallel Simulation Engine for Molecular Dynamics,"
16 * J. Comput. Chem. 26, pp. 252-271 (2005))
17 *
18 * 2. Redistributions of source code must retain the above copyright
19 * notice, this list of conditions and the following disclaimer.
20 *
21 * 3. Redistributions in binary form must reproduce the above copyright
22 * notice, this list of conditions and the following disclaimer in the
23 * documentation and/or other materials provided with the
24 * distribution.
25 *
26 * This software is provided "AS IS," without a warranty of any
27 * kind. All express or implied conditions, representations and
28 * warranties, including any implied warranty of merchantability,
29 * fitness for a particular purpose or non-infringement, are hereby
30 * excluded. The University of Notre Dame and its licensors shall not
31 * be liable for any damages suffered by licensee as a result of
32 * using, modifying or distributing the software or its
33 * derivatives. In no event will the University of Notre Dame or its
34 * licensors be liable for any lost revenue, profit or data, or for
35 * direct, indirect, special, consequential, incidental or punitive
36 * damages, however caused and regardless of the theory of liability,
37 * arising out of the use of or inability to use software, even if the
38 * University of Notre Dame has been advised of the possibility of
39 * such damages.
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/fElectrostaticSummationMethod.h"
59 #include "UseTheForce/DarkSide/fElectrostaticScreeningMethod.h"
60 #include "UseTheForce/DarkSide/fSwitchingFunctionType.h"
61 #include "UseTheForce/doForces_interface.h"
62 #include "UseTheForce/DarkSide/neighborLists_interface.h"
63 #include "UseTheForce/DarkSide/electrostatic_interface.h"
64 #include "UseTheForce/DarkSide/switcheroo_interface.h"
65 #include "utils/MemoryUtils.hpp"
66 #include "utils/simError.h"
67 #include "selection/SelectionManager.hpp"
68 #include "io/ForceFieldOptions.hpp"
69 #include "UseTheForce/ForceField.hpp"
70
71
72 #ifdef IS_MPI
73 #include "UseTheForce/mpiComponentPlan.h"
74 #include "UseTheForce/DarkSide/simParallel_interface.h"
75 #endif
76
77 namespace oopse {
78 std::set<int> getRigidSet(int index, std::map<int, std::set<int> >& container) {
79 std::map<int, std::set<int> >::iterator i = container.find(index);
80 std::set<int> result;
81 if (i != container.end()) {
82 result = i->second;
83 }
84
85 return result;
86 }
87
88 SimInfo::SimInfo(ForceField* ff, Globals* simParams) :
89 forceField_(ff), simParams_(simParams),
90 ndf_(0), fdf_local(0), ndfRaw_(0), ndfTrans_(0), nZconstraint_(0),
91 nGlobalMols_(0), nGlobalAtoms_(0), nGlobalCutoffGroups_(0),
92 nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0),
93 nAtoms_(0), nBonds_(0), nBends_(0), nTorsions_(0), nInversions_(0),
94 nRigidBodies_(0), nIntegrableObjects_(0), nCutoffGroups_(0),
95 nConstraints_(0), sman_(NULL), fortranInitialized_(false),
96 calcBoxDipole_(false), useAtomicVirial_(true) {
97
98
99 MoleculeStamp* molStamp;
100 int nMolWithSameStamp;
101 int nCutoffAtoms = 0; // number of atoms belong to cutoff groups
102 int nGroups = 0; //total cutoff groups defined in meta-data file
103 CutoffGroupStamp* cgStamp;
104 RigidBodyStamp* rbStamp;
105 int nRigidAtoms = 0;
106
107 std::vector<Component*> components = simParams->getComponents();
108
109 for (std::vector<Component*>::iterator i = components.begin(); i !=components.end(); ++i) {
110 molStamp = (*i)->getMoleculeStamp();
111 nMolWithSameStamp = (*i)->getNMol();
112
113 addMoleculeStamp(molStamp, nMolWithSameStamp);
114
115 //calculate atoms in molecules
116 nGlobalAtoms_ += molStamp->getNAtoms() *nMolWithSameStamp;
117
118 //calculate atoms in cutoff groups
119 int nAtomsInGroups = 0;
120 int nCutoffGroupsInStamp = molStamp->getNCutoffGroups();
121
122 for (int j=0; j < nCutoffGroupsInStamp; j++) {
123 cgStamp = molStamp->getCutoffGroupStamp(j);
124 nAtomsInGroups += cgStamp->getNMembers();
125 }
126
127 nGroups += nCutoffGroupsInStamp * nMolWithSameStamp;
128
129 nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp;
130
131 //calculate atoms in rigid bodies
132 int nAtomsInRigidBodies = 0;
133 int nRigidBodiesInStamp = molStamp->getNRigidBodies();
134
135 for (int j=0; j < nRigidBodiesInStamp; j++) {
136 rbStamp = molStamp->getRigidBodyStamp(j);
137 nAtomsInRigidBodies += rbStamp->getNMembers();
138 }
139
140 nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp;
141 nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp;
142
143 }
144
145 //every free atom (atom does not belong to cutoff groups) is a cutoff
146 //group therefore the total number of cutoff groups in the system is
147 //equal to the total number of atoms minus number of atoms belong to
148 //cutoff group defined in meta-data file plus the number of cutoff
149 //groups defined in meta-data file
150 nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
151
152 //every free atom (atom does not belong to rigid bodies) is an
153 //integrable object therefore the total number of integrable objects
154 //in the system is equal to the total number of atoms minus number of
155 //atoms belong to rigid body defined in meta-data file plus the number
156 //of rigid bodies defined in meta-data file
157 nGlobalIntegrableObjects_ = nGlobalAtoms_ - nRigidAtoms
158 + nGlobalRigidBodies_;
159
160 nGlobalMols_ = molStampIds_.size();
161 molToProcMap_.resize(nGlobalMols_);
162 }
163
164 SimInfo::~SimInfo() {
165 std::map<int, Molecule*>::iterator i;
166 for (i = molecules_.begin(); i != molecules_.end(); ++i) {
167 delete i->second;
168 }
169 molecules_.clear();
170
171 delete sman_;
172 delete simParams_;
173 delete forceField_;
174 }
175
176 int SimInfo::getNGlobalConstraints() {
177 int nGlobalConstraints;
178 #ifdef IS_MPI
179 MPI_Allreduce(&nConstraints_, &nGlobalConstraints, 1, MPI_INT, MPI_SUM,
180 MPI_COMM_WORLD);
181 #else
182 nGlobalConstraints = nConstraints_;
183 #endif
184 return nGlobalConstraints;
185 }
186
187 bool SimInfo::addMolecule(Molecule* mol) {
188 MoleculeIterator i;
189
190 i = molecules_.find(mol->getGlobalIndex());
191 if (i == molecules_.end() ) {
192
193 molecules_.insert(std::make_pair(mol->getGlobalIndex(), mol));
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 addInteractionPairs(mol);
206
207 return true;
208 } else {
209 return false;
210 }
211 }
212
213 bool SimInfo::removeMolecule(Molecule* mol) {
214 MoleculeIterator i;
215 i = molecules_.find(mol->getGlobalIndex());
216
217 if (i != molecules_.end() ) {
218
219 assert(mol == i->second);
220
221 nAtoms_ -= mol->getNAtoms();
222 nBonds_ -= mol->getNBonds();
223 nBends_ -= mol->getNBends();
224 nTorsions_ -= mol->getNTorsions();
225 nInversions_ -= mol->getNInversions();
226 nRigidBodies_ -= mol->getNRigidBodies();
227 nIntegrableObjects_ -= mol->getNIntegrableObjects();
228 nCutoffGroups_ -= mol->getNCutoffGroups();
229 nConstraints_ -= mol->getNConstraintPairs();
230
231 removeInteractionPairs(mol);
232 molecules_.erase(mol->getGlobalIndex());
233
234 delete mol;
235
236 return true;
237 } else {
238 return false;
239 }
240
241
242 }
243
244
245 Molecule* SimInfo::beginMolecule(MoleculeIterator& i) {
246 i = molecules_.begin();
247 return i == molecules_.end() ? NULL : i->second;
248 }
249
250 Molecule* SimInfo::nextMolecule(MoleculeIterator& i) {
251 ++i;
252 return i == molecules_.end() ? NULL : i->second;
253 }
254
255
256 void SimInfo::calcNdf() {
257 int ndf_local;
258 MoleculeIterator i;
259 std::vector<StuntDouble*>::iterator j;
260 Molecule* mol;
261 StuntDouble* integrableObject;
262
263 ndf_local = 0;
264
265 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
266 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
267 integrableObject = mol->nextIntegrableObject(j)) {
268
269 ndf_local += 3;
270
271 if (integrableObject->isDirectional()) {
272 if (integrableObject->isLinear()) {
273 ndf_local += 2;
274 } else {
275 ndf_local += 3;
276 }
277 }
278
279 }
280 }
281
282 // n_constraints is local, so subtract them on each processor
283 ndf_local -= nConstraints_;
284
285 #ifdef IS_MPI
286 MPI_Allreduce(&ndf_local,&ndf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
287 #else
288 ndf_ = ndf_local;
289 #endif
290
291 // nZconstraints_ is global, as are the 3 COM translations for the
292 // entire system:
293 ndf_ = ndf_ - 3 - nZconstraint_;
294
295 }
296
297 int SimInfo::getFdf() {
298 #ifdef IS_MPI
299 MPI_Allreduce(&fdf_local,&fdf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
300 #else
301 fdf_ = fdf_local;
302 #endif
303 return fdf_;
304 }
305
306 void SimInfo::calcNdfRaw() {
307 int ndfRaw_local;
308
309 MoleculeIterator i;
310 std::vector<StuntDouble*>::iterator j;
311 Molecule* mol;
312 StuntDouble* integrableObject;
313
314 // Raw degrees of freedom that we have to set
315 ndfRaw_local = 0;
316
317 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
318 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
319 integrableObject = mol->nextIntegrableObject(j)) {
320
321 ndfRaw_local += 3;
322
323 if (integrableObject->isDirectional()) {
324 if (integrableObject->isLinear()) {
325 ndfRaw_local += 2;
326 } else {
327 ndfRaw_local += 3;
328 }
329 }
330
331 }
332 }
333
334 #ifdef IS_MPI
335 MPI_Allreduce(&ndfRaw_local,&ndfRaw_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
336 #else
337 ndfRaw_ = ndfRaw_local;
338 #endif
339 }
340
341 void SimInfo::calcNdfTrans() {
342 int ndfTrans_local;
343
344 ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_;
345
346
347 #ifdef IS_MPI
348 MPI_Allreduce(&ndfTrans_local,&ndfTrans_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
349 #else
350 ndfTrans_ = ndfTrans_local;
351 #endif
352
353 ndfTrans_ = ndfTrans_ - 3 - nZconstraint_;
354
355 }
356
357 void SimInfo::addInteractionPairs(Molecule* mol) {
358 ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
359 std::vector<Bond*>::iterator bondIter;
360 std::vector<Bend*>::iterator bendIter;
361 std::vector<Torsion*>::iterator torsionIter;
362 std::vector<Inversion*>::iterator inversionIter;
363 Bond* bond;
364 Bend* bend;
365 Torsion* torsion;
366 Inversion* inversion;
367 int a;
368 int b;
369 int c;
370 int d;
371
372 // atomGroups can be used to add special interaction maps between
373 // groups of atoms that are in two separate rigid bodies.
374 // However, most site-site interactions between two rigid bodies
375 // are probably not special, just the ones between the physically
376 // bonded atoms. Interactions *within* a single rigid body should
377 // always be excluded. These are done at the bottom of this
378 // function.
379
380 std::map<int, std::set<int> > atomGroups;
381 Molecule::RigidBodyIterator rbIter;
382 RigidBody* rb;
383 Molecule::IntegrableObjectIterator ii;
384 StuntDouble* integrableObject;
385
386 for (integrableObject = mol->beginIntegrableObject(ii);
387 integrableObject != NULL;
388 integrableObject = mol->nextIntegrableObject(ii)) {
389
390 if (integrableObject->isRigidBody()) {
391 rb = static_cast<RigidBody*>(integrableObject);
392 std::vector<Atom*> atoms = rb->getAtoms();
393 std::set<int> rigidAtoms;
394 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
395 rigidAtoms.insert(atoms[i]->getGlobalIndex());
396 }
397 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
398 atomGroups.insert(std::map<int, std::set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
399 }
400 } else {
401 std::set<int> oneAtomSet;
402 oneAtomSet.insert(integrableObject->getGlobalIndex());
403 atomGroups.insert(std::map<int, std::set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));
404 }
405 }
406
407 for (bond= mol->beginBond(bondIter); bond != NULL;
408 bond = mol->nextBond(bondIter)) {
409
410 a = bond->getAtomA()->getGlobalIndex();
411 b = bond->getAtomB()->getGlobalIndex();
412
413 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
414 oneTwoInteractions_.addPair(a, b);
415 } else {
416 excludedInteractions_.addPair(a, b);
417 }
418 }
419
420 for (bend= mol->beginBend(bendIter); bend != NULL;
421 bend = mol->nextBend(bendIter)) {
422
423 a = bend->getAtomA()->getGlobalIndex();
424 b = bend->getAtomB()->getGlobalIndex();
425 c = bend->getAtomC()->getGlobalIndex();
426
427 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
428 oneTwoInteractions_.addPair(a, b);
429 oneTwoInteractions_.addPair(b, c);
430 } else {
431 excludedInteractions_.addPair(a, b);
432 excludedInteractions_.addPair(b, c);
433 }
434
435 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
436 oneThreeInteractions_.addPair(a, c);
437 } else {
438 excludedInteractions_.addPair(a, c);
439 }
440 }
441
442 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL;
443 torsion = mol->nextTorsion(torsionIter)) {
444
445 a = torsion->getAtomA()->getGlobalIndex();
446 b = torsion->getAtomB()->getGlobalIndex();
447 c = torsion->getAtomC()->getGlobalIndex();
448 d = torsion->getAtomD()->getGlobalIndex();
449
450 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
451 oneTwoInteractions_.addPair(a, b);
452 oneTwoInteractions_.addPair(b, c);
453 oneTwoInteractions_.addPair(c, d);
454 } else {
455 excludedInteractions_.addPair(a, b);
456 excludedInteractions_.addPair(b, c);
457 excludedInteractions_.addPair(c, d);
458 }
459
460 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
461 oneThreeInteractions_.addPair(a, c);
462 oneThreeInteractions_.addPair(b, d);
463 } else {
464 excludedInteractions_.addPair(a, c);
465 excludedInteractions_.addPair(b, d);
466 }
467
468 if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
469 oneFourInteractions_.addPair(a, d);
470 } else {
471 excludedInteractions_.addPair(a, d);
472 }
473 }
474
475 for (inversion= mol->beginInversion(inversionIter); inversion != NULL;
476 inversion = mol->nextInversion(inversionIter)) {
477
478 a = inversion->getAtomA()->getGlobalIndex();
479 b = inversion->getAtomB()->getGlobalIndex();
480 c = inversion->getAtomC()->getGlobalIndex();
481 d = inversion->getAtomD()->getGlobalIndex();
482
483 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
484 oneTwoInteractions_.addPair(a, b);
485 oneTwoInteractions_.addPair(a, c);
486 oneTwoInteractions_.addPair(a, d);
487 } else {
488 excludedInteractions_.addPair(a, b);
489 excludedInteractions_.addPair(a, c);
490 excludedInteractions_.addPair(a, d);
491 }
492
493 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
494 oneThreeInteractions_.addPair(b, c);
495 oneThreeInteractions_.addPair(b, d);
496 oneThreeInteractions_.addPair(c, d);
497 } else {
498 excludedInteractions_.addPair(b, c);
499 excludedInteractions_.addPair(b, d);
500 excludedInteractions_.addPair(c, d);
501 }
502 }
503
504 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
505 rb = mol->nextRigidBody(rbIter)) {
506 std::vector<Atom*> atoms = rb->getAtoms();
507 for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) {
508 for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
509 a = atoms[i]->getGlobalIndex();
510 b = atoms[j]->getGlobalIndex();
511 excludedInteractions_.addPair(a, b);
512 }
513 }
514 }
515
516 }
517
518 void SimInfo::removeInteractionPairs(Molecule* mol) {
519 ForceFieldOptions& options_ = forceField_->getForceFieldOptions();
520 std::vector<Bond*>::iterator bondIter;
521 std::vector<Bend*>::iterator bendIter;
522 std::vector<Torsion*>::iterator torsionIter;
523 std::vector<Inversion*>::iterator inversionIter;
524 Bond* bond;
525 Bend* bend;
526 Torsion* torsion;
527 Inversion* inversion;
528 int a;
529 int b;
530 int c;
531 int d;
532
533 std::map<int, std::set<int> > atomGroups;
534 Molecule::RigidBodyIterator rbIter;
535 RigidBody* rb;
536 Molecule::IntegrableObjectIterator ii;
537 StuntDouble* integrableObject;
538
539 for (integrableObject = mol->beginIntegrableObject(ii);
540 integrableObject != NULL;
541 integrableObject = mol->nextIntegrableObject(ii)) {
542
543 if (integrableObject->isRigidBody()) {
544 rb = static_cast<RigidBody*>(integrableObject);
545 std::vector<Atom*> atoms = rb->getAtoms();
546 std::set<int> rigidAtoms;
547 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
548 rigidAtoms.insert(atoms[i]->getGlobalIndex());
549 }
550 for (int i = 0; i < static_cast<int>(atoms.size()); ++i) {
551 atomGroups.insert(std::map<int, std::set<int> >::value_type(atoms[i]->getGlobalIndex(), rigidAtoms));
552 }
553 } else {
554 std::set<int> oneAtomSet;
555 oneAtomSet.insert(integrableObject->getGlobalIndex());
556 atomGroups.insert(std::map<int, std::set<int> >::value_type(integrableObject->getGlobalIndex(), oneAtomSet));
557 }
558 }
559
560 for (bond= mol->beginBond(bondIter); bond != NULL;
561 bond = mol->nextBond(bondIter)) {
562
563 a = bond->getAtomA()->getGlobalIndex();
564 b = bond->getAtomB()->getGlobalIndex();
565
566 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
567 oneTwoInteractions_.removePair(a, b);
568 } else {
569 excludedInteractions_.removePair(a, b);
570 }
571 }
572
573 for (bend= mol->beginBend(bendIter); bend != NULL;
574 bend = mol->nextBend(bendIter)) {
575
576 a = bend->getAtomA()->getGlobalIndex();
577 b = bend->getAtomB()->getGlobalIndex();
578 c = bend->getAtomC()->getGlobalIndex();
579
580 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
581 oneTwoInteractions_.removePair(a, b);
582 oneTwoInteractions_.removePair(b, c);
583 } else {
584 excludedInteractions_.removePair(a, b);
585 excludedInteractions_.removePair(b, c);
586 }
587
588 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
589 oneThreeInteractions_.removePair(a, c);
590 } else {
591 excludedInteractions_.removePair(a, c);
592 }
593 }
594
595 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL;
596 torsion = mol->nextTorsion(torsionIter)) {
597
598 a = torsion->getAtomA()->getGlobalIndex();
599 b = torsion->getAtomB()->getGlobalIndex();
600 c = torsion->getAtomC()->getGlobalIndex();
601 d = torsion->getAtomD()->getGlobalIndex();
602
603 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
604 oneTwoInteractions_.removePair(a, b);
605 oneTwoInteractions_.removePair(b, c);
606 oneTwoInteractions_.removePair(c, d);
607 } else {
608 excludedInteractions_.removePair(a, b);
609 excludedInteractions_.removePair(b, c);
610 excludedInteractions_.removePair(c, d);
611 }
612
613 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
614 oneThreeInteractions_.removePair(a, c);
615 oneThreeInteractions_.removePair(b, d);
616 } else {
617 excludedInteractions_.removePair(a, c);
618 excludedInteractions_.removePair(b, d);
619 }
620
621 if (options_.havevdw14scale() || options_.haveelectrostatic14scale()) {
622 oneFourInteractions_.removePair(a, d);
623 } else {
624 excludedInteractions_.removePair(a, d);
625 }
626 }
627
628 for (inversion= mol->beginInversion(inversionIter); inversion != NULL;
629 inversion = mol->nextInversion(inversionIter)) {
630
631 a = inversion->getAtomA()->getGlobalIndex();
632 b = inversion->getAtomB()->getGlobalIndex();
633 c = inversion->getAtomC()->getGlobalIndex();
634 d = inversion->getAtomD()->getGlobalIndex();
635
636 if (options_.havevdw12scale() || options_.haveelectrostatic12scale()) {
637 oneTwoInteractions_.removePair(a, b);
638 oneTwoInteractions_.removePair(a, c);
639 oneTwoInteractions_.removePair(a, d);
640 } else {
641 excludedInteractions_.removePair(a, b);
642 excludedInteractions_.removePair(a, c);
643 excludedInteractions_.removePair(a, d);
644 }
645
646 if (options_.havevdw13scale() || options_.haveelectrostatic13scale()) {
647 oneThreeInteractions_.removePair(b, c);
648 oneThreeInteractions_.removePair(b, d);
649 oneThreeInteractions_.removePair(c, d);
650 } else {
651 excludedInteractions_.removePair(b, c);
652 excludedInteractions_.removePair(b, d);
653 excludedInteractions_.removePair(c, d);
654 }
655 }
656
657 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
658 rb = mol->nextRigidBody(rbIter)) {
659 std::vector<Atom*> atoms = rb->getAtoms();
660 for (int i = 0; i < static_cast<int>(atoms.size()) -1 ; ++i) {
661 for (int j = i + 1; j < static_cast<int>(atoms.size()); ++j) {
662 a = atoms[i]->getGlobalIndex();
663 b = atoms[j]->getGlobalIndex();
664 excludedInteractions_.removePair(a, b);
665 }
666 }
667 }
668
669 }
670
671
672 void SimInfo::addMoleculeStamp(MoleculeStamp* molStamp, int nmol) {
673 int curStampId;
674
675 //index from 0
676 curStampId = moleculeStamps_.size();
677
678 moleculeStamps_.push_back(molStamp);
679 molStampIds_.insert(molStampIds_.end(), nmol, curStampId);
680 }
681
682 void SimInfo::update() {
683
684 setupSimType();
685
686 #ifdef IS_MPI
687 setupFortranParallel();
688 #endif
689
690 setupFortranSim();
691
692 //setup fortran force field
693 /** @deprecate */
694 int isError = 0;
695
696 setupCutoff();
697
698 setupElectrostaticSummationMethod( isError );
699 setupSwitchingFunction();
700 setupAccumulateBoxDipole();
701
702 if(isError){
703 sprintf( painCave.errMsg,
704 "ForceField error: There was an error initializing the forceField in fortran.\n" );
705 painCave.isFatal = 1;
706 simError();
707 }
708
709 calcNdf();
710 calcNdfRaw();
711 calcNdfTrans();
712
713 fortranInitialized_ = true;
714 }
715
716 std::set<AtomType*> SimInfo::getUniqueAtomTypes() {
717 SimInfo::MoleculeIterator mi;
718 Molecule* mol;
719 Molecule::AtomIterator ai;
720 Atom* atom;
721 std::set<AtomType*> atomTypes;
722
723 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
724
725 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
726 atomTypes.insert(atom->getAtomType());
727 }
728
729 }
730
731 return atomTypes;
732 }
733
734 void SimInfo::setupSimType() {
735 std::set<AtomType*>::iterator i;
736 std::set<AtomType*> atomTypes;
737 atomTypes = getUniqueAtomTypes();
738
739 int useLennardJones = 0;
740 int useElectrostatic = 0;
741 int useEAM = 0;
742 int useSC = 0;
743 int useCharge = 0;
744 int useDirectional = 0;
745 int useDipole = 0;
746 int useGayBerne = 0;
747 int useSticky = 0;
748 int useStickyPower = 0;
749 int useShape = 0;
750 int useFLARB = 0; //it is not in AtomType yet
751 int useDirectionalAtom = 0;
752 int useElectrostatics = 0;
753 //usePBC and useRF are from simParams
754 int usePBC = simParams_->getUsePeriodicBoundaryConditions();
755 int useRF;
756 int useSF;
757 int useSP;
758 int useBoxDipole;
759
760 std::string myMethod;
761
762 // set the useRF logical
763 useRF = 0;
764 useSF = 0;
765 useSP = 0;
766
767
768 if (simParams_->haveElectrostaticSummationMethod()) {
769 std::string myMethod = simParams_->getElectrostaticSummationMethod();
770 toUpper(myMethod);
771 if (myMethod == "REACTION_FIELD"){
772 useRF = 1;
773 } else if (myMethod == "SHIFTED_FORCE"){
774 useSF = 1;
775 } else if (myMethod == "SHIFTED_POTENTIAL"){
776 useSP = 1;
777 }
778 }
779
780 if (simParams_->haveAccumulateBoxDipole())
781 if (simParams_->getAccumulateBoxDipole())
782 useBoxDipole = 1;
783
784 useAtomicVirial_ = simParams_->getUseAtomicVirial();
785
786 //loop over all of the atom types
787 for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
788 useLennardJones |= (*i)->isLennardJones();
789 useElectrostatic |= (*i)->isElectrostatic();
790 useEAM |= (*i)->isEAM();
791 useSC |= (*i)->isSC();
792 useCharge |= (*i)->isCharge();
793 useDirectional |= (*i)->isDirectional();
794 useDipole |= (*i)->isDipole();
795 useGayBerne |= (*i)->isGayBerne();
796 useSticky |= (*i)->isSticky();
797 useStickyPower |= (*i)->isStickyPower();
798 useShape |= (*i)->isShape();
799 }
800
801 if (useSticky || useStickyPower || useDipole || useGayBerne || useShape) {
802 useDirectionalAtom = 1;
803 }
804
805 if (useCharge || useDipole) {
806 useElectrostatics = 1;
807 }
808
809 #ifdef IS_MPI
810 int temp;
811
812 temp = usePBC;
813 MPI_Allreduce(&temp, &usePBC, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
814
815 temp = useDirectionalAtom;
816 MPI_Allreduce(&temp, &useDirectionalAtom, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
817
818 temp = useLennardJones;
819 MPI_Allreduce(&temp, &useLennardJones, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
820
821 temp = useElectrostatics;
822 MPI_Allreduce(&temp, &useElectrostatics, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
823
824 temp = useCharge;
825 MPI_Allreduce(&temp, &useCharge, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
826
827 temp = useDipole;
828 MPI_Allreduce(&temp, &useDipole, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
829
830 temp = useSticky;
831 MPI_Allreduce(&temp, &useSticky, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
832
833 temp = useStickyPower;
834 MPI_Allreduce(&temp, &useStickyPower, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
835
836 temp = useGayBerne;
837 MPI_Allreduce(&temp, &useGayBerne, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
838
839 temp = useEAM;
840 MPI_Allreduce(&temp, &useEAM, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
841
842 temp = useSC;
843 MPI_Allreduce(&temp, &useSC, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
844
845 temp = useShape;
846 MPI_Allreduce(&temp, &useShape, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
847
848 temp = useFLARB;
849 MPI_Allreduce(&temp, &useFLARB, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
850
851 temp = useRF;
852 MPI_Allreduce(&temp, &useRF, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
853
854 temp = useSF;
855 MPI_Allreduce(&temp, &useSF, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
856
857 temp = useSP;
858 MPI_Allreduce(&temp, &useSP, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
859
860 temp = useBoxDipole;
861 MPI_Allreduce(&temp, &useBoxDipole, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
862
863 temp = useAtomicVirial_;
864 MPI_Allreduce(&temp, &useAtomicVirial_, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
865
866 #endif
867
868 fInfo_.SIM_uses_PBC = usePBC;
869 fInfo_.SIM_uses_DirectionalAtoms = useDirectionalAtom;
870 fInfo_.SIM_uses_LennardJones = useLennardJones;
871 fInfo_.SIM_uses_Electrostatics = useElectrostatics;
872 fInfo_.SIM_uses_Charges = useCharge;
873 fInfo_.SIM_uses_Dipoles = useDipole;
874 fInfo_.SIM_uses_Sticky = useSticky;
875 fInfo_.SIM_uses_StickyPower = useStickyPower;
876 fInfo_.SIM_uses_GayBerne = useGayBerne;
877 fInfo_.SIM_uses_EAM = useEAM;
878 fInfo_.SIM_uses_SC = useSC;
879 fInfo_.SIM_uses_Shapes = useShape;
880 fInfo_.SIM_uses_FLARB = useFLARB;
881 fInfo_.SIM_uses_RF = useRF;
882 fInfo_.SIM_uses_SF = useSF;
883 fInfo_.SIM_uses_SP = useSP;
884 fInfo_.SIM_uses_BoxDipole = useBoxDipole;
885 fInfo_.SIM_uses_AtomicVirial = useAtomicVirial_;
886 }
887
888 void SimInfo::setupFortranSim() {
889 int isError;
890 int nExclude, nOneTwo, nOneThree, nOneFour;
891 std::vector<int> fortranGlobalGroupMembership;
892
893 isError = 0;
894
895 //globalGroupMembership_ is filled by SimCreator
896 for (int i = 0; i < nGlobalAtoms_; i++) {
897 fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1);
898 }
899
900 //calculate mass ratio of cutoff group
901 std::vector<RealType> mfact;
902 SimInfo::MoleculeIterator mi;
903 Molecule* mol;
904 Molecule::CutoffGroupIterator ci;
905 CutoffGroup* cg;
906 Molecule::AtomIterator ai;
907 Atom* atom;
908 RealType totalMass;
909
910 //to avoid memory reallocation, reserve enough space for mfact
911 mfact.reserve(getNCutoffGroups());
912
913 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
914 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
915
916 totalMass = cg->getMass();
917 for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
918 // Check for massless groups - set mfact to 1 if true
919 if (totalMass != 0)
920 mfact.push_back(atom->getMass()/totalMass);
921 else
922 mfact.push_back( 1.0 );
923 }
924 }
925 }
926
927 //fill ident array of local atoms (it is actually ident of AtomType, it is so confusing !!!)
928 std::vector<int> identArray;
929
930 //to avoid memory reallocation, reserve enough space identArray
931 identArray.reserve(getNAtoms());
932
933 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
934 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
935 identArray.push_back(atom->getIdent());
936 }
937 }
938
939 //fill molMembershipArray
940 //molMembershipArray is filled by SimCreator
941 std::vector<int> molMembershipArray(nGlobalAtoms_);
942 for (int i = 0; i < nGlobalAtoms_; i++) {
943 molMembershipArray[i] = globalMolMembership_[i] + 1;
944 }
945
946 //setup fortran simulation
947
948 nExclude = excludedInteractions_.getSize();
949 nOneTwo = oneTwoInteractions_.getSize();
950 nOneThree = oneThreeInteractions_.getSize();
951 nOneFour = oneFourInteractions_.getSize();
952
953 std::cerr << "exculdes:\n";
954 std::cerr << excludedInteractions_;
955 std::cerr << "\noneTwo:\n";
956 std::cerr << oneTwoInteractions_;
957 std::cerr << "\noneThree:\n";
958 std::cerr << oneThreeInteractions_;
959 std::cerr << "\noneFour:\n";
960 std::cerr << oneFourInteractions_;
961
962 int* excludeList = excludedInteractions_.getPairList();
963 int* oneTwoList = oneTwoInteractions_.getPairList();
964 int* oneThreeList = oneThreeInteractions_.getPairList();
965 int* oneFourList = oneFourInteractions_.getPairList();
966
967 setFortranSim( &fInfo_, &nGlobalAtoms_, &nAtoms_, &identArray[0],
968 &nExclude, excludeList,
969 &nOneTwo, oneTwoList,
970 &nOneThree, oneThreeList,
971 &nOneFour, oneFourList,
972 &molMembershipArray[0], &mfact[0], &nCutoffGroups_,
973 &fortranGlobalGroupMembership[0], &isError);
974
975 if( isError ){
976
977 sprintf( painCave.errMsg,
978 "There was an error setting the simulation information in fortran.\n" );
979 painCave.isFatal = 1;
980 painCave.severity = OOPSE_ERROR;
981 simError();
982 }
983
984
985 sprintf( checkPointMsg,
986 "succesfully sent the simulation information to fortran.\n");
987
988 errorCheckPoint();
989
990 // Setup number of neighbors in neighbor list if present
991 if (simParams_->haveNeighborListNeighbors()) {
992 int nlistNeighbors = simParams_->getNeighborListNeighbors();
993 setNeighbors(&nlistNeighbors);
994 }
995
996
997 }
998
999
1000 void SimInfo::setupFortranParallel() {
1001 #ifdef IS_MPI
1002 //SimInfo is responsible for creating localToGlobalAtomIndex and localToGlobalGroupIndex
1003 std::vector<int> localToGlobalAtomIndex(getNAtoms(), 0);
1004 std::vector<int> localToGlobalCutoffGroupIndex;
1005 SimInfo::MoleculeIterator mi;
1006 Molecule::AtomIterator ai;
1007 Molecule::CutoffGroupIterator ci;
1008 Molecule* mol;
1009 Atom* atom;
1010 CutoffGroup* cg;
1011 mpiSimData parallelData;
1012 int isError;
1013
1014 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
1015
1016 //local index(index in DataStorge) of atom is important
1017 for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
1018 localToGlobalAtomIndex[atom->getLocalIndex()] = atom->getGlobalIndex() + 1;
1019 }
1020
1021 //local index of cutoff group is trivial, it only depends on the order of travesing
1022 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
1023 localToGlobalCutoffGroupIndex.push_back(cg->getGlobalIndex() + 1);
1024 }
1025
1026 }
1027
1028 //fill up mpiSimData struct
1029 parallelData.nMolGlobal = getNGlobalMolecules();
1030 parallelData.nMolLocal = getNMolecules();
1031 parallelData.nAtomsGlobal = getNGlobalAtoms();
1032 parallelData.nAtomsLocal = getNAtoms();
1033 parallelData.nGroupsGlobal = getNGlobalCutoffGroups();
1034 parallelData.nGroupsLocal = getNCutoffGroups();
1035 parallelData.myNode = worldRank;
1036 MPI_Comm_size(MPI_COMM_WORLD, &(parallelData.nProcessors));
1037
1038 //pass mpiSimData struct and index arrays to fortran
1039 setFsimParallel(&parallelData, &(parallelData.nAtomsLocal),
1040 &localToGlobalAtomIndex[0], &(parallelData.nGroupsLocal),
1041 &localToGlobalCutoffGroupIndex[0], &isError);
1042
1043 if (isError) {
1044 sprintf(painCave.errMsg,
1045 "mpiRefresh errror: fortran didn't like something we gave it.\n");
1046 painCave.isFatal = 1;
1047 simError();
1048 }
1049
1050 sprintf(checkPointMsg, " mpiRefresh successful.\n");
1051 errorCheckPoint();
1052
1053 #endif
1054 }
1055
1056 void SimInfo::setupCutoff() {
1057
1058 ForceFieldOptions& forceFieldOptions_ = forceField_->getForceFieldOptions();
1059
1060 // Check the cutoff policy
1061 int cp = TRADITIONAL_CUTOFF_POLICY; // Set to traditional by default
1062
1063 // Set LJ shifting bools to false
1064 ljsp_ = false;
1065 ljsf_ = false;
1066
1067 std::string myPolicy;
1068 if (forceFieldOptions_.haveCutoffPolicy()){
1069 myPolicy = forceFieldOptions_.getCutoffPolicy();
1070 }else if (simParams_->haveCutoffPolicy()) {
1071 myPolicy = simParams_->getCutoffPolicy();
1072 }
1073
1074 if (!myPolicy.empty()){
1075 toUpper(myPolicy);
1076 if (myPolicy == "MIX") {
1077 cp = MIX_CUTOFF_POLICY;
1078 } else {
1079 if (myPolicy == "MAX") {
1080 cp = MAX_CUTOFF_POLICY;
1081 } else {
1082 if (myPolicy == "TRADITIONAL") {
1083 cp = TRADITIONAL_CUTOFF_POLICY;
1084 } else {
1085 // throw error
1086 sprintf( painCave.errMsg,
1087 "SimInfo error: Unknown cutoffPolicy. (Input file specified %s .)\n\tcutoffPolicy must be one of: \"Mix\", \"Max\", or \"Traditional\".", myPolicy.c_str() );
1088 painCave.isFatal = 1;
1089 simError();
1090 }
1091 }
1092 }
1093 }
1094 notifyFortranCutoffPolicy(&cp);
1095
1096 // Check the Skin Thickness for neighborlists
1097 RealType skin;
1098 if (simParams_->haveSkinThickness()) {
1099 skin = simParams_->getSkinThickness();
1100 notifyFortranSkinThickness(&skin);
1101 }
1102
1103 // Check if the cutoff was set explicitly:
1104 if (simParams_->haveCutoffRadius()) {
1105 rcut_ = simParams_->getCutoffRadius();
1106 if (simParams_->haveSwitchingRadius()) {
1107 rsw_ = simParams_->getSwitchingRadius();
1108 } else {
1109 if (fInfo_.SIM_uses_Charges |
1110 fInfo_.SIM_uses_Dipoles |
1111 fInfo_.SIM_uses_RF) {
1112
1113 rsw_ = 0.85 * rcut_;
1114 sprintf(painCave.errMsg,
1115 "SimCreator Warning: No value was set for the switchingRadius.\n"
1116 "\tOOPSE will use a default value of 85 percent of the cutoffRadius.\n"
1117 "\tswitchingRadius = %f. for this simulation\n", rsw_);
1118 painCave.isFatal = 0;
1119 simError();
1120 } else {
1121 rsw_ = rcut_;
1122 sprintf(painCave.errMsg,
1123 "SimCreator Warning: No value was set for the switchingRadius.\n"
1124 "\tOOPSE will use the same value as the cutoffRadius.\n"
1125 "\tswitchingRadius = %f. for this simulation\n", rsw_);
1126 painCave.isFatal = 0;
1127 simError();
1128 }
1129 }
1130
1131 if (simParams_->haveElectrostaticSummationMethod()) {
1132 std::string myMethod = simParams_->getElectrostaticSummationMethod();
1133 toUpper(myMethod);
1134
1135 if (myMethod == "SHIFTED_POTENTIAL") {
1136 ljsp_ = true;
1137 } else if (myMethod == "SHIFTED_FORCE") {
1138 ljsf_ = true;
1139 }
1140 }
1141 notifyFortranCutoffs(&rcut_, &rsw_, &ljsp_, &ljsf_);
1142
1143 } else {
1144
1145 // For electrostatic atoms, we'll assume a large safe value:
1146 if (fInfo_.SIM_uses_Charges | fInfo_.SIM_uses_Dipoles | fInfo_.SIM_uses_RF) {
1147 sprintf(painCave.errMsg,
1148 "SimCreator Warning: No value was set for the cutoffRadius.\n"
1149 "\tOOPSE will use a default value of 15.0 angstroms"
1150 "\tfor the cutoffRadius.\n");
1151 painCave.isFatal = 0;
1152 simError();
1153 rcut_ = 15.0;
1154
1155 if (simParams_->haveElectrostaticSummationMethod()) {
1156 std::string myMethod = simParams_->getElectrostaticSummationMethod();
1157 toUpper(myMethod);
1158
1159 // For the time being, we're tethering the LJ shifted behavior to the
1160 // electrostaticSummationMethod keyword options
1161 if (myMethod == "SHIFTED_POTENTIAL") {
1162 ljsp_ = true;
1163 } else if (myMethod == "SHIFTED_FORCE") {
1164 ljsf_ = true;
1165 }
1166 if (myMethod == "SHIFTED_POTENTIAL" || myMethod == "SHIFTED_FORCE") {
1167 if (simParams_->haveSwitchingRadius()){
1168 sprintf(painCave.errMsg,
1169 "SimInfo Warning: A value was set for the switchingRadius\n"
1170 "\teven though the electrostaticSummationMethod was\n"
1171 "\tset to %s\n", myMethod.c_str());
1172 painCave.isFatal = 1;
1173 simError();
1174 }
1175 }
1176 }
1177
1178 if (simParams_->haveSwitchingRadius()){
1179 rsw_ = simParams_->getSwitchingRadius();
1180 } else {
1181 sprintf(painCave.errMsg,
1182 "SimCreator Warning: No value was set for switchingRadius.\n"
1183 "\tOOPSE will use a default value of\n"
1184 "\t0.85 * cutoffRadius for the switchingRadius\n");
1185 painCave.isFatal = 0;
1186 simError();
1187 rsw_ = 0.85 * rcut_;
1188 }
1189
1190 notifyFortranCutoffs(&rcut_, &rsw_, &ljsp_, &ljsf_);
1191
1192 } else {
1193 // We didn't set rcut explicitly, and we don't have electrostatic atoms, so
1194 // We'll punt and let fortran figure out the cutoffs later.
1195
1196 notifyFortranYouAreOnYourOwn();
1197
1198 }
1199 }
1200 }
1201
1202 void SimInfo::setupElectrostaticSummationMethod( int isError ) {
1203
1204 int errorOut;
1205 int esm = NONE;
1206 int sm = UNDAMPED;
1207 RealType alphaVal;
1208 RealType dielectric;
1209
1210 errorOut = isError;
1211
1212 if (simParams_->haveElectrostaticSummationMethod()) {
1213 std::string myMethod = simParams_->getElectrostaticSummationMethod();
1214 toUpper(myMethod);
1215 if (myMethod == "NONE") {
1216 esm = NONE;
1217 } else {
1218 if (myMethod == "SWITCHING_FUNCTION") {
1219 esm = SWITCHING_FUNCTION;
1220 } else {
1221 if (myMethod == "SHIFTED_POTENTIAL") {
1222 esm = SHIFTED_POTENTIAL;
1223 } else {
1224 if (myMethod == "SHIFTED_FORCE") {
1225 esm = SHIFTED_FORCE;
1226 } else {
1227 if (myMethod == "REACTION_FIELD") {
1228 esm = REACTION_FIELD;
1229 dielectric = simParams_->getDielectric();
1230 if (!simParams_->haveDielectric()) {
1231 // throw warning
1232 sprintf( painCave.errMsg,
1233 "SimInfo warning: dielectric was not specified in the input file\n\tfor the reaction field correction method.\n"
1234 "\tA default value of %f will be used for the dielectric.\n", dielectric);
1235 painCave.isFatal = 0;
1236 simError();
1237 }
1238 } else {
1239 // throw error
1240 sprintf( painCave.errMsg,
1241 "SimInfo error: Unknown electrostaticSummationMethod.\n"
1242 "\t(Input file specified %s .)\n"
1243 "\telectrostaticSummationMethod must be one of: \"none\",\n"
1244 "\t\"shifted_potential\", \"shifted_force\", or \n"
1245 "\t\"reaction_field\".\n", myMethod.c_str() );
1246 painCave.isFatal = 1;
1247 simError();
1248 }
1249 }
1250 }
1251 }
1252 }
1253 }
1254
1255 if (simParams_->haveElectrostaticScreeningMethod()) {
1256 std::string myScreen = simParams_->getElectrostaticScreeningMethod();
1257 toUpper(myScreen);
1258 if (myScreen == "UNDAMPED") {
1259 sm = UNDAMPED;
1260 } else {
1261 if (myScreen == "DAMPED") {
1262 sm = DAMPED;
1263 if (!simParams_->haveDampingAlpha()) {
1264 // first set a cutoff dependent alpha value
1265 // we assume alpha depends linearly with rcut from 0 to 20.5 ang
1266 alphaVal = 0.5125 - rcut_* 0.025;
1267 // for values rcut > 20.5, alpha is zero
1268 if (alphaVal < 0) alphaVal = 0;
1269
1270 // throw warning
1271 sprintf( painCave.errMsg,
1272 "SimInfo warning: dampingAlpha was not specified in the input file.\n"
1273 "\tA default value of %f (1/ang) will be used for the cutoff of\n\t%f (ang).\n", alphaVal, rcut_);
1274 painCave.isFatal = 0;
1275 simError();
1276 } else {
1277 alphaVal = simParams_->getDampingAlpha();
1278 }
1279
1280 } else {
1281 // throw error
1282 sprintf( painCave.errMsg,
1283 "SimInfo error: Unknown electrostaticScreeningMethod.\n"
1284 "\t(Input file specified %s .)\n"
1285 "\telectrostaticScreeningMethod must be one of: \"undamped\"\n"
1286 "or \"damped\".\n", myScreen.c_str() );
1287 painCave.isFatal = 1;
1288 simError();
1289 }
1290 }
1291 }
1292
1293 // let's pass some summation method variables to fortran
1294 setElectrostaticSummationMethod( &esm );
1295 setFortranElectrostaticMethod( &esm );
1296 setScreeningMethod( &sm );
1297 setDampingAlpha( &alphaVal );
1298 setReactionFieldDielectric( &dielectric );
1299 initFortranFF( &errorOut );
1300 }
1301
1302 void SimInfo::setupSwitchingFunction() {
1303 int ft = CUBIC;
1304
1305 if (simParams_->haveSwitchingFunctionType()) {
1306 std::string funcType = simParams_->getSwitchingFunctionType();
1307 toUpper(funcType);
1308 if (funcType == "CUBIC") {
1309 ft = CUBIC;
1310 } else {
1311 if (funcType == "FIFTH_ORDER_POLYNOMIAL") {
1312 ft = FIFTH_ORDER_POLY;
1313 } else {
1314 // throw error
1315 sprintf( painCave.errMsg,
1316 "SimInfo error: Unknown switchingFunctionType. (Input file specified %s .)\n\tswitchingFunctionType must be one of: \"cubic\" or \"fifth_order_polynomial\".", funcType.c_str() );
1317 painCave.isFatal = 1;
1318 simError();
1319 }
1320 }
1321 }
1322
1323 // send switching function notification to switcheroo
1324 setFunctionType(&ft);
1325
1326 }
1327
1328 void SimInfo::setupAccumulateBoxDipole() {
1329
1330 // we only call setAccumulateBoxDipole if the accumulateBoxDipole parameter is true
1331 if ( simParams_->haveAccumulateBoxDipole() )
1332 if ( simParams_->getAccumulateBoxDipole() ) {
1333 setAccumulateBoxDipole();
1334 calcBoxDipole_ = true;
1335 }
1336
1337 }
1338
1339 void SimInfo::addProperty(GenericData* genData) {
1340 properties_.addProperty(genData);
1341 }
1342
1343 void SimInfo::removeProperty(const std::string& propName) {
1344 properties_.removeProperty(propName);
1345 }
1346
1347 void SimInfo::clearProperties() {
1348 properties_.clearProperties();
1349 }
1350
1351 std::vector<std::string> SimInfo::getPropertyNames() {
1352 return properties_.getPropertyNames();
1353 }
1354
1355 std::vector<GenericData*> SimInfo::getProperties() {
1356 return properties_.getProperties();
1357 }
1358
1359 GenericData* SimInfo::getPropertyByName(const std::string& propName) {
1360 return properties_.getPropertyByName(propName);
1361 }
1362
1363 void SimInfo::setSnapshotManager(SnapshotManager* sman) {
1364 if (sman_ == sman) {
1365 return;
1366 }
1367 delete sman_;
1368 sman_ = sman;
1369
1370 Molecule* mol;
1371 RigidBody* rb;
1372 Atom* atom;
1373 SimInfo::MoleculeIterator mi;
1374 Molecule::RigidBodyIterator rbIter;
1375 Molecule::AtomIterator atomIter;;
1376
1377 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
1378
1379 for (atom = mol->beginAtom(atomIter); atom != NULL; atom = mol->nextAtom(atomIter)) {
1380 atom->setSnapshotManager(sman_);
1381 }
1382
1383 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
1384 rb->setSnapshotManager(sman_);
1385 }
1386 }
1387
1388 }
1389
1390 Vector3d SimInfo::getComVel(){
1391 SimInfo::MoleculeIterator i;
1392 Molecule* mol;
1393
1394 Vector3d comVel(0.0);
1395 RealType totalMass = 0.0;
1396
1397
1398 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1399 RealType mass = mol->getMass();
1400 totalMass += mass;
1401 comVel += mass * mol->getComVel();
1402 }
1403
1404 #ifdef IS_MPI
1405 RealType tmpMass = totalMass;
1406 Vector3d tmpComVel(comVel);
1407 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1408 MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1409 #endif
1410
1411 comVel /= totalMass;
1412
1413 return comVel;
1414 }
1415
1416 Vector3d SimInfo::getCom(){
1417 SimInfo::MoleculeIterator i;
1418 Molecule* mol;
1419
1420 Vector3d com(0.0);
1421 RealType totalMass = 0.0;
1422
1423 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1424 RealType mass = mol->getMass();
1425 totalMass += mass;
1426 com += mass * mol->getCom();
1427 }
1428
1429 #ifdef IS_MPI
1430 RealType tmpMass = totalMass;
1431 Vector3d tmpCom(com);
1432 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1433 MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1434 #endif
1435
1436 com /= totalMass;
1437
1438 return com;
1439
1440 }
1441
1442 std::ostream& operator <<(std::ostream& o, SimInfo& info) {
1443
1444 return o;
1445 }
1446
1447
1448 /*
1449 Returns center of mass and center of mass velocity in one function call.
1450 */
1451
1452 void SimInfo::getComAll(Vector3d &com, Vector3d &comVel){
1453 SimInfo::MoleculeIterator i;
1454 Molecule* mol;
1455
1456
1457 RealType totalMass = 0.0;
1458
1459
1460 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1461 RealType mass = mol->getMass();
1462 totalMass += mass;
1463 com += mass * mol->getCom();
1464 comVel += mass * mol->getComVel();
1465 }
1466
1467 #ifdef IS_MPI
1468 RealType tmpMass = totalMass;
1469 Vector3d tmpCom(com);
1470 Vector3d tmpComVel(comVel);
1471 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1472 MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1473 MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1474 #endif
1475
1476 com /= totalMass;
1477 comVel /= totalMass;
1478 }
1479
1480 /*
1481 Return intertia tensor for entire system and angular momentum Vector.
1482
1483
1484 [ Ixx -Ixy -Ixz ]
1485 J =| -Iyx Iyy -Iyz |
1486 [ -Izx -Iyz Izz ]
1487 */
1488
1489 void SimInfo::getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum){
1490
1491
1492 RealType xx = 0.0;
1493 RealType yy = 0.0;
1494 RealType zz = 0.0;
1495 RealType xy = 0.0;
1496 RealType xz = 0.0;
1497 RealType yz = 0.0;
1498 Vector3d com(0.0);
1499 Vector3d comVel(0.0);
1500
1501 getComAll(com, comVel);
1502
1503 SimInfo::MoleculeIterator i;
1504 Molecule* mol;
1505
1506 Vector3d thisq(0.0);
1507 Vector3d thisv(0.0);
1508
1509 RealType thisMass = 0.0;
1510
1511
1512
1513
1514 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1515
1516 thisq = mol->getCom()-com;
1517 thisv = mol->getComVel()-comVel;
1518 thisMass = mol->getMass();
1519 // Compute moment of intertia coefficients.
1520 xx += thisq[0]*thisq[0]*thisMass;
1521 yy += thisq[1]*thisq[1]*thisMass;
1522 zz += thisq[2]*thisq[2]*thisMass;
1523
1524 // compute products of intertia
1525 xy += thisq[0]*thisq[1]*thisMass;
1526 xz += thisq[0]*thisq[2]*thisMass;
1527 yz += thisq[1]*thisq[2]*thisMass;
1528
1529 angularMomentum += cross( thisq, thisv ) * thisMass;
1530
1531 }
1532
1533
1534 inertiaTensor(0,0) = yy + zz;
1535 inertiaTensor(0,1) = -xy;
1536 inertiaTensor(0,2) = -xz;
1537 inertiaTensor(1,0) = -xy;
1538 inertiaTensor(1,1) = xx + zz;
1539 inertiaTensor(1,2) = -yz;
1540 inertiaTensor(2,0) = -xz;
1541 inertiaTensor(2,1) = -yz;
1542 inertiaTensor(2,2) = xx + yy;
1543
1544 #ifdef IS_MPI
1545 Mat3x3d tmpI(inertiaTensor);
1546 Vector3d tmpAngMom;
1547 MPI_Allreduce(tmpI.getArrayPointer(), inertiaTensor.getArrayPointer(),9,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1548 MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1549 #endif
1550
1551 return;
1552 }
1553
1554 //Returns the angular momentum of the system
1555 Vector3d SimInfo::getAngularMomentum(){
1556
1557 Vector3d com(0.0);
1558 Vector3d comVel(0.0);
1559 Vector3d angularMomentum(0.0);
1560
1561 getComAll(com,comVel);
1562
1563 SimInfo::MoleculeIterator i;
1564 Molecule* mol;
1565
1566 Vector3d thisr(0.0);
1567 Vector3d thisp(0.0);
1568
1569 RealType thisMass;
1570
1571 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1572 thisMass = mol->getMass();
1573 thisr = mol->getCom()-com;
1574 thisp = (mol->getComVel()-comVel)*thisMass;
1575
1576 angularMomentum += cross( thisr, thisp );
1577
1578 }
1579
1580 #ifdef IS_MPI
1581 Vector3d tmpAngMom;
1582 MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_REALTYPE,MPI_SUM, MPI_COMM_WORLD);
1583 #endif
1584
1585 return angularMomentum;
1586 }
1587
1588 StuntDouble* SimInfo::getIOIndexToIntegrableObject(int index) {
1589 return IOIndexToIntegrableObject.at(index);
1590 }
1591
1592 void SimInfo::setIOIndexToIntegrableObject(const std::vector<StuntDouble*>& v) {
1593 IOIndexToIntegrableObject= v;
1594 }
1595
1596 /* Returns the Volume of the simulation based on a ellipsoid with semi-axes
1597 based on the radius of gyration V=4/3*Pi*R_1*R_2*R_3
1598 where R_i are related to the principle inertia moments R_i = sqrt(C*I_i/N), this reduces to
1599 V = 4/3*Pi*(C/N)^3/2*sqrt(det(I)). See S.E. Baltazar et. al. Comp. Mat. Sci. 37 (2006) 526-536.
1600 */
1601 void SimInfo::getGyrationalVolume(RealType &volume){
1602 Mat3x3d intTensor;
1603 RealType det;
1604 Vector3d dummyAngMom;
1605 RealType sysconstants;
1606 RealType geomCnst;
1607
1608 geomCnst = 3.0/2.0;
1609 /* Get the inertial tensor and angular momentum for free*/
1610 getInertiaTensor(intTensor,dummyAngMom);
1611
1612 det = intTensor.determinant();
1613 sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_;
1614 volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(det);
1615 return;
1616 }
1617
1618 void SimInfo::getGyrationalVolume(RealType &volume, RealType &detI){
1619 Mat3x3d intTensor;
1620 Vector3d dummyAngMom;
1621 RealType sysconstants;
1622 RealType geomCnst;
1623
1624 geomCnst = 3.0/2.0;
1625 /* Get the inertial tensor and angular momentum for free*/
1626 getInertiaTensor(intTensor,dummyAngMom);
1627
1628 detI = intTensor.determinant();
1629 sysconstants = geomCnst/(RealType)nGlobalIntegrableObjects_;
1630 volume = 4.0/3.0*NumericConstant::PI*pow(sysconstants,3.0/2.0)*sqrt(detI);
1631 return;
1632 }
1633 /*
1634 void SimInfo::setStuntDoubleFromGlobalIndex(std::vector<StuntDouble*> v) {
1635 assert( v.size() == nAtoms_ + nRigidBodies_);
1636 sdByGlobalIndex_ = v;
1637 }
1638
1639 StuntDouble* SimInfo::getStuntDoubleFromGlobalIndex(int index) {
1640 //assert(index < nAtoms_ + nRigidBodies_);
1641 return sdByGlobalIndex_.at(index);
1642 }
1643 */
1644 }//end namespace oopse
1645