ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/brains/SimInfo.cpp
Revision: 2309
Committed: Sun Sep 18 20:45:38 2005 UTC (18 years, 10 months ago) by chrisfen
File size: 34919 byte(s)
Log Message:
reaction field seems to work now, still need to do some testing...

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
52 #include "brains/SimInfo.hpp"
53 #include "math/Vector3.hpp"
54 #include "primitives/Molecule.hpp"
55 #include "UseTheForce/fCutoffPolicy.h"
56 #include "UseTheForce/DarkSide/fElectrostaticSummationMethod.h"
57 #include "UseTheForce/doForces_interface.h"
58 #include "UseTheForce/DarkSide/electrostatic_interface.h"
59 #include "UseTheForce/notifyCutoffs_interface.h"
60 #include "utils/MemoryUtils.hpp"
61 #include "utils/simError.h"
62 #include "selection/SelectionManager.hpp"
63
64 #ifdef IS_MPI
65 #include "UseTheForce/mpiComponentPlan.h"
66 #include "UseTheForce/DarkSide/simParallel_interface.h"
67 #endif
68
69 namespace oopse {
70
71 SimInfo::SimInfo(MakeStamps* stamps, std::vector<std::pair<MoleculeStamp*, int> >& molStampPairs,
72 ForceField* ff, Globals* simParams) :
73 stamps_(stamps), forceField_(ff), simParams_(simParams),
74 ndf_(0), ndfRaw_(0), ndfTrans_(0), nZconstraint_(0),
75 nGlobalMols_(0), nGlobalAtoms_(0), nGlobalCutoffGroups_(0),
76 nGlobalIntegrableObjects_(0), nGlobalRigidBodies_(0),
77 nAtoms_(0), nBonds_(0), nBends_(0), nTorsions_(0), nRigidBodies_(0),
78 nIntegrableObjects_(0), nCutoffGroups_(0), nConstraints_(0),
79 sman_(NULL), fortranInitialized_(false) {
80
81
82 std::vector<std::pair<MoleculeStamp*, int> >::iterator i;
83 MoleculeStamp* molStamp;
84 int nMolWithSameStamp;
85 int nCutoffAtoms = 0; // number of atoms belong to cutoff groups
86 int nGroups = 0; //total cutoff groups defined in meta-data file
87 CutoffGroupStamp* cgStamp;
88 RigidBodyStamp* rbStamp;
89 int nRigidAtoms = 0;
90
91 for (i = molStampPairs.begin(); i !=molStampPairs.end(); ++i) {
92 molStamp = i->first;
93 nMolWithSameStamp = i->second;
94
95 addMoleculeStamp(molStamp, nMolWithSameStamp);
96
97 //calculate atoms in molecules
98 nGlobalAtoms_ += molStamp->getNAtoms() *nMolWithSameStamp;
99
100
101 //calculate atoms in cutoff groups
102 int nAtomsInGroups = 0;
103 int nCutoffGroupsInStamp = molStamp->getNCutoffGroups();
104
105 for (int j=0; j < nCutoffGroupsInStamp; j++) {
106 cgStamp = molStamp->getCutoffGroup(j);
107 nAtomsInGroups += cgStamp->getNMembers();
108 }
109
110 nGroups += nCutoffGroupsInStamp * nMolWithSameStamp;
111 nCutoffAtoms += nAtomsInGroups * nMolWithSameStamp;
112
113 //calculate atoms in rigid bodies
114 int nAtomsInRigidBodies = 0;
115 int nRigidBodiesInStamp = molStamp->getNRigidBodies();
116
117 for (int j=0; j < nRigidBodiesInStamp; j++) {
118 rbStamp = molStamp->getRigidBody(j);
119 nAtomsInRigidBodies += rbStamp->getNMembers();
120 }
121
122 nGlobalRigidBodies_ += nRigidBodiesInStamp * nMolWithSameStamp;
123 nRigidAtoms += nAtomsInRigidBodies * nMolWithSameStamp;
124
125 }
126
127 //every free atom (atom does not belong to cutoff groups) is a cutoff group
128 //therefore the total number of cutoff groups in the system is equal to
129 //the total number of atoms minus number of atoms belong to cutoff group defined in meta-data
130 //file plus the number of cutoff groups defined in meta-data file
131 nGlobalCutoffGroups_ = nGlobalAtoms_ - nCutoffAtoms + nGroups;
132
133 //every free atom (atom does not belong to rigid bodies) is an integrable object
134 //therefore the total number of integrable objects in the system is equal to
135 //the total number of atoms minus number of atoms belong to rigid body defined in meta-data
136 //file plus the number of rigid bodies defined in meta-data file
137 nGlobalIntegrableObjects_ = nGlobalAtoms_ - nRigidAtoms + nGlobalRigidBodies_;
138
139 nGlobalMols_ = molStampIds_.size();
140
141 #ifdef IS_MPI
142 molToProcMap_.resize(nGlobalMols_);
143 #endif
144
145 }
146
147 SimInfo::~SimInfo() {
148 std::map<int, Molecule*>::iterator i;
149 for (i = molecules_.begin(); i != molecules_.end(); ++i) {
150 delete i->second;
151 }
152 molecules_.clear();
153
154 delete stamps_;
155 delete sman_;
156 delete simParams_;
157 delete forceField_;
158 }
159
160 int SimInfo::getNGlobalConstraints() {
161 int nGlobalConstraints;
162 #ifdef IS_MPI
163 MPI_Allreduce(&nConstraints_, &nGlobalConstraints, 1, MPI_INT, MPI_SUM,
164 MPI_COMM_WORLD);
165 #else
166 nGlobalConstraints = nConstraints_;
167 #endif
168 return nGlobalConstraints;
169 }
170
171 bool SimInfo::addMolecule(Molecule* mol) {
172 MoleculeIterator i;
173
174 i = molecules_.find(mol->getGlobalIndex());
175 if (i == molecules_.end() ) {
176
177 molecules_.insert(std::make_pair(mol->getGlobalIndex(), mol));
178
179 nAtoms_ += mol->getNAtoms();
180 nBonds_ += mol->getNBonds();
181 nBends_ += mol->getNBends();
182 nTorsions_ += mol->getNTorsions();
183 nRigidBodies_ += mol->getNRigidBodies();
184 nIntegrableObjects_ += mol->getNIntegrableObjects();
185 nCutoffGroups_ += mol->getNCutoffGroups();
186 nConstraints_ += mol->getNConstraintPairs();
187
188 addExcludePairs(mol);
189
190 return true;
191 } else {
192 return false;
193 }
194 }
195
196 bool SimInfo::removeMolecule(Molecule* mol) {
197 MoleculeIterator i;
198 i = molecules_.find(mol->getGlobalIndex());
199
200 if (i != molecules_.end() ) {
201
202 assert(mol == i->second);
203
204 nAtoms_ -= mol->getNAtoms();
205 nBonds_ -= mol->getNBonds();
206 nBends_ -= mol->getNBends();
207 nTorsions_ -= mol->getNTorsions();
208 nRigidBodies_ -= mol->getNRigidBodies();
209 nIntegrableObjects_ -= mol->getNIntegrableObjects();
210 nCutoffGroups_ -= mol->getNCutoffGroups();
211 nConstraints_ -= mol->getNConstraintPairs();
212
213 removeExcludePairs(mol);
214 molecules_.erase(mol->getGlobalIndex());
215
216 delete mol;
217
218 return true;
219 } else {
220 return false;
221 }
222
223
224 }
225
226
227 Molecule* SimInfo::beginMolecule(MoleculeIterator& i) {
228 i = molecules_.begin();
229 return i == molecules_.end() ? NULL : i->second;
230 }
231
232 Molecule* SimInfo::nextMolecule(MoleculeIterator& i) {
233 ++i;
234 return i == molecules_.end() ? NULL : i->second;
235 }
236
237
238 void SimInfo::calcNdf() {
239 int ndf_local;
240 MoleculeIterator i;
241 std::vector<StuntDouble*>::iterator j;
242 Molecule* mol;
243 StuntDouble* integrableObject;
244
245 ndf_local = 0;
246
247 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
248 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
249 integrableObject = mol->nextIntegrableObject(j)) {
250
251 ndf_local += 3;
252
253 if (integrableObject->isDirectional()) {
254 if (integrableObject->isLinear()) {
255 ndf_local += 2;
256 } else {
257 ndf_local += 3;
258 }
259 }
260
261 }//end for (integrableObject)
262 }// end for (mol)
263
264 // n_constraints is local, so subtract them on each processor
265 ndf_local -= nConstraints_;
266
267 #ifdef IS_MPI
268 MPI_Allreduce(&ndf_local,&ndf_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
269 #else
270 ndf_ = ndf_local;
271 #endif
272
273 // nZconstraints_ is global, as are the 3 COM translations for the
274 // entire system:
275 ndf_ = ndf_ - 3 - nZconstraint_;
276
277 }
278
279 void SimInfo::calcNdfRaw() {
280 int ndfRaw_local;
281
282 MoleculeIterator i;
283 std::vector<StuntDouble*>::iterator j;
284 Molecule* mol;
285 StuntDouble* integrableObject;
286
287 // Raw degrees of freedom that we have to set
288 ndfRaw_local = 0;
289
290 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
291 for (integrableObject = mol->beginIntegrableObject(j); integrableObject != NULL;
292 integrableObject = mol->nextIntegrableObject(j)) {
293
294 ndfRaw_local += 3;
295
296 if (integrableObject->isDirectional()) {
297 if (integrableObject->isLinear()) {
298 ndfRaw_local += 2;
299 } else {
300 ndfRaw_local += 3;
301 }
302 }
303
304 }
305 }
306
307 #ifdef IS_MPI
308 MPI_Allreduce(&ndfRaw_local,&ndfRaw_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
309 #else
310 ndfRaw_ = ndfRaw_local;
311 #endif
312 }
313
314 void SimInfo::calcNdfTrans() {
315 int ndfTrans_local;
316
317 ndfTrans_local = 3 * nIntegrableObjects_ - nConstraints_;
318
319
320 #ifdef IS_MPI
321 MPI_Allreduce(&ndfTrans_local,&ndfTrans_,1,MPI_INT,MPI_SUM, MPI_COMM_WORLD);
322 #else
323 ndfTrans_ = ndfTrans_local;
324 #endif
325
326 ndfTrans_ = ndfTrans_ - 3 - nZconstraint_;
327
328 }
329
330 void SimInfo::addExcludePairs(Molecule* mol) {
331 std::vector<Bond*>::iterator bondIter;
332 std::vector<Bend*>::iterator bendIter;
333 std::vector<Torsion*>::iterator torsionIter;
334 Bond* bond;
335 Bend* bend;
336 Torsion* torsion;
337 int a;
338 int b;
339 int c;
340 int d;
341
342 for (bond= mol->beginBond(bondIter); bond != NULL; bond = mol->nextBond(bondIter)) {
343 a = bond->getAtomA()->getGlobalIndex();
344 b = bond->getAtomB()->getGlobalIndex();
345 exclude_.addPair(a, b);
346 }
347
348 for (bend= mol->beginBend(bendIter); bend != NULL; bend = mol->nextBend(bendIter)) {
349 a = bend->getAtomA()->getGlobalIndex();
350 b = bend->getAtomB()->getGlobalIndex();
351 c = bend->getAtomC()->getGlobalIndex();
352
353 exclude_.addPair(a, b);
354 exclude_.addPair(a, c);
355 exclude_.addPair(b, c);
356 }
357
358 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL; torsion = mol->nextTorsion(torsionIter)) {
359 a = torsion->getAtomA()->getGlobalIndex();
360 b = torsion->getAtomB()->getGlobalIndex();
361 c = torsion->getAtomC()->getGlobalIndex();
362 d = torsion->getAtomD()->getGlobalIndex();
363
364 exclude_.addPair(a, b);
365 exclude_.addPair(a, c);
366 exclude_.addPair(a, d);
367 exclude_.addPair(b, c);
368 exclude_.addPair(b, d);
369 exclude_.addPair(c, d);
370 }
371
372 Molecule::RigidBodyIterator rbIter;
373 RigidBody* rb;
374 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
375 std::vector<Atom*> atoms = rb->getAtoms();
376 for (int i = 0; i < atoms.size() -1 ; ++i) {
377 for (int j = i + 1; j < atoms.size(); ++j) {
378 a = atoms[i]->getGlobalIndex();
379 b = atoms[j]->getGlobalIndex();
380 exclude_.addPair(a, b);
381 }
382 }
383 }
384
385 }
386
387 void SimInfo::removeExcludePairs(Molecule* mol) {
388 std::vector<Bond*>::iterator bondIter;
389 std::vector<Bend*>::iterator bendIter;
390 std::vector<Torsion*>::iterator torsionIter;
391 Bond* bond;
392 Bend* bend;
393 Torsion* torsion;
394 int a;
395 int b;
396 int c;
397 int d;
398
399 for (bond= mol->beginBond(bondIter); bond != NULL; bond = mol->nextBond(bondIter)) {
400 a = bond->getAtomA()->getGlobalIndex();
401 b = bond->getAtomB()->getGlobalIndex();
402 exclude_.removePair(a, b);
403 }
404
405 for (bend= mol->beginBend(bendIter); bend != NULL; bend = mol->nextBend(bendIter)) {
406 a = bend->getAtomA()->getGlobalIndex();
407 b = bend->getAtomB()->getGlobalIndex();
408 c = bend->getAtomC()->getGlobalIndex();
409
410 exclude_.removePair(a, b);
411 exclude_.removePair(a, c);
412 exclude_.removePair(b, c);
413 }
414
415 for (torsion= mol->beginTorsion(torsionIter); torsion != NULL; torsion = mol->nextTorsion(torsionIter)) {
416 a = torsion->getAtomA()->getGlobalIndex();
417 b = torsion->getAtomB()->getGlobalIndex();
418 c = torsion->getAtomC()->getGlobalIndex();
419 d = torsion->getAtomD()->getGlobalIndex();
420
421 exclude_.removePair(a, b);
422 exclude_.removePair(a, c);
423 exclude_.removePair(a, d);
424 exclude_.removePair(b, c);
425 exclude_.removePair(b, d);
426 exclude_.removePair(c, d);
427 }
428
429 Molecule::RigidBodyIterator rbIter;
430 RigidBody* rb;
431 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
432 std::vector<Atom*> atoms = rb->getAtoms();
433 for (int i = 0; i < atoms.size() -1 ; ++i) {
434 for (int j = i + 1; j < atoms.size(); ++j) {
435 a = atoms[i]->getGlobalIndex();
436 b = atoms[j]->getGlobalIndex();
437 exclude_.removePair(a, b);
438 }
439 }
440 }
441
442 }
443
444
445 void SimInfo::addMoleculeStamp(MoleculeStamp* molStamp, int nmol) {
446 int curStampId;
447
448 //index from 0
449 curStampId = moleculeStamps_.size();
450
451 moleculeStamps_.push_back(molStamp);
452 molStampIds_.insert(molStampIds_.end(), nmol, curStampId);
453 }
454
455 void SimInfo::update() {
456
457 setupSimType();
458
459 #ifdef IS_MPI
460 setupFortranParallel();
461 #endif
462
463 setupFortranSim();
464
465 //setup fortran force field
466 /** @deprecate */
467 int isError = 0;
468
469 setupElectrostaticSummationMethod( isError );
470
471 if(isError){
472 sprintf( painCave.errMsg,
473 "ForceField error: There was an error initializing the forceField in fortran.\n" );
474 painCave.isFatal = 1;
475 simError();
476 }
477
478
479 setupCutoff();
480
481 calcNdf();
482 calcNdfRaw();
483 calcNdfTrans();
484
485 fortranInitialized_ = true;
486 }
487
488 std::set<AtomType*> SimInfo::getUniqueAtomTypes() {
489 SimInfo::MoleculeIterator mi;
490 Molecule* mol;
491 Molecule::AtomIterator ai;
492 Atom* atom;
493 std::set<AtomType*> atomTypes;
494
495 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
496
497 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
498 atomTypes.insert(atom->getAtomType());
499 }
500
501 }
502
503 return atomTypes;
504 }
505
506 void SimInfo::setupSimType() {
507 std::set<AtomType*>::iterator i;
508 std::set<AtomType*> atomTypes;
509 atomTypes = getUniqueAtomTypes();
510
511 int useLennardJones = 0;
512 int useElectrostatic = 0;
513 int useEAM = 0;
514 int useCharge = 0;
515 int useDirectional = 0;
516 int useDipole = 0;
517 int useGayBerne = 0;
518 int useSticky = 0;
519 int useStickyPower = 0;
520 int useShape = 0;
521 int useFLARB = 0; //it is not in AtomType yet
522 int useDirectionalAtom = 0;
523 int useElectrostatics = 0;
524 //usePBC and useRF are from simParams
525 int usePBC = simParams_->getPBC();
526
527 //loop over all of the atom types
528 for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
529 useLennardJones |= (*i)->isLennardJones();
530 useElectrostatic |= (*i)->isElectrostatic();
531 useEAM |= (*i)->isEAM();
532 useCharge |= (*i)->isCharge();
533 useDirectional |= (*i)->isDirectional();
534 useDipole |= (*i)->isDipole();
535 useGayBerne |= (*i)->isGayBerne();
536 useSticky |= (*i)->isSticky();
537 useStickyPower |= (*i)->isStickyPower();
538 useShape |= (*i)->isShape();
539 }
540
541 if (useSticky || useStickyPower || useDipole || useGayBerne || useShape) {
542 useDirectionalAtom = 1;
543 }
544
545 if (useCharge || useDipole) {
546 useElectrostatics = 1;
547 }
548
549 #ifdef IS_MPI
550 int temp;
551
552 temp = usePBC;
553 MPI_Allreduce(&temp, &usePBC, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
554
555 temp = useDirectionalAtom;
556 MPI_Allreduce(&temp, &useDirectionalAtom, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
557
558 temp = useLennardJones;
559 MPI_Allreduce(&temp, &useLennardJones, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
560
561 temp = useElectrostatics;
562 MPI_Allreduce(&temp, &useElectrostatics, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
563
564 temp = useCharge;
565 MPI_Allreduce(&temp, &useCharge, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
566
567 temp = useDipole;
568 MPI_Allreduce(&temp, &useDipole, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
569
570 temp = useSticky;
571 MPI_Allreduce(&temp, &useSticky, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
572
573 temp = useStickyPower;
574 MPI_Allreduce(&temp, &useStickyPower, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
575
576 temp = useGayBerne;
577 MPI_Allreduce(&temp, &useGayBerne, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
578
579 temp = useEAM;
580 MPI_Allreduce(&temp, &useEAM, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
581
582 temp = useShape;
583 MPI_Allreduce(&temp, &useShape, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
584
585 temp = useFLARB;
586 MPI_Allreduce(&temp, &useFLARB, 1, MPI_INT, MPI_LOR, MPI_COMM_WORLD);
587
588 #endif
589
590 fInfo_.SIM_uses_PBC = usePBC;
591 fInfo_.SIM_uses_DirectionalAtoms = useDirectionalAtom;
592 fInfo_.SIM_uses_LennardJones = useLennardJones;
593 fInfo_.SIM_uses_Electrostatics = useElectrostatics;
594 fInfo_.SIM_uses_Charges = useCharge;
595 fInfo_.SIM_uses_Dipoles = useDipole;
596 fInfo_.SIM_uses_Sticky = useSticky;
597 fInfo_.SIM_uses_StickyPower = useStickyPower;
598 fInfo_.SIM_uses_GayBerne = useGayBerne;
599 fInfo_.SIM_uses_EAM = useEAM;
600 fInfo_.SIM_uses_Shapes = useShape;
601 fInfo_.SIM_uses_FLARB = useFLARB;
602
603 if( fInfo_.SIM_uses_Dipoles && fInfo_.SIM_uses_RF) {
604
605 if (simParams_->haveDielectric()) {
606 fInfo_.dielect = simParams_->getDielectric();
607 } else {
608 sprintf(painCave.errMsg,
609 "SimSetup Error: No Dielectric constant was set.\n"
610 "\tYou are trying to use Reaction Field without"
611 "\tsetting a dielectric constant!\n");
612 painCave.isFatal = 1;
613 simError();
614 }
615
616 } else {
617 fInfo_.dielect = 0.0;
618 }
619
620 }
621
622 void SimInfo::setupFortranSim() {
623 int isError;
624 int nExclude;
625 std::vector<int> fortranGlobalGroupMembership;
626
627 nExclude = exclude_.getSize();
628 isError = 0;
629
630 //globalGroupMembership_ is filled by SimCreator
631 for (int i = 0; i < nGlobalAtoms_; i++) {
632 fortranGlobalGroupMembership.push_back(globalGroupMembership_[i] + 1);
633 }
634
635 //calculate mass ratio of cutoff group
636 std::vector<double> mfact;
637 SimInfo::MoleculeIterator mi;
638 Molecule* mol;
639 Molecule::CutoffGroupIterator ci;
640 CutoffGroup* cg;
641 Molecule::AtomIterator ai;
642 Atom* atom;
643 double totalMass;
644
645 //to avoid memory reallocation, reserve enough space for mfact
646 mfact.reserve(getNCutoffGroups());
647
648 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
649 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
650
651 totalMass = cg->getMass();
652 for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
653 mfact.push_back(atom->getMass()/totalMass);
654 }
655
656 }
657 }
658
659 //fill ident array of local atoms (it is actually ident of AtomType, it is so confusing !!!)
660 std::vector<int> identArray;
661
662 //to avoid memory reallocation, reserve enough space identArray
663 identArray.reserve(getNAtoms());
664
665 for(mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
666 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
667 identArray.push_back(atom->getIdent());
668 }
669 }
670
671 //fill molMembershipArray
672 //molMembershipArray is filled by SimCreator
673 std::vector<int> molMembershipArray(nGlobalAtoms_);
674 for (int i = 0; i < nGlobalAtoms_; i++) {
675 molMembershipArray[i] = globalMolMembership_[i] + 1;
676 }
677
678 //setup fortran simulation
679 int nGlobalExcludes = 0;
680 int* globalExcludes = NULL;
681 int* excludeList = exclude_.getExcludeList();
682 setFortranSim( &fInfo_, &nGlobalAtoms_, &nAtoms_, &identArray[0], &nExclude, excludeList ,
683 &nGlobalExcludes, globalExcludes, &molMembershipArray[0],
684 &mfact[0], &nCutoffGroups_, &fortranGlobalGroupMembership[0], &isError);
685
686 if( isError ){
687
688 sprintf( painCave.errMsg,
689 "There was an error setting the simulation information in fortran.\n" );
690 painCave.isFatal = 1;
691 painCave.severity = OOPSE_ERROR;
692 simError();
693 }
694
695 #ifdef IS_MPI
696 sprintf( checkPointMsg,
697 "succesfully sent the simulation information to fortran.\n");
698 MPIcheckPoint();
699 #endif // is_mpi
700 }
701
702
703 #ifdef IS_MPI
704 void SimInfo::setupFortranParallel() {
705
706 //SimInfo is responsible for creating localToGlobalAtomIndex and localToGlobalGroupIndex
707 std::vector<int> localToGlobalAtomIndex(getNAtoms(), 0);
708 std::vector<int> localToGlobalCutoffGroupIndex;
709 SimInfo::MoleculeIterator mi;
710 Molecule::AtomIterator ai;
711 Molecule::CutoffGroupIterator ci;
712 Molecule* mol;
713 Atom* atom;
714 CutoffGroup* cg;
715 mpiSimData parallelData;
716 int isError;
717
718 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
719
720 //local index(index in DataStorge) of atom is important
721 for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
722 localToGlobalAtomIndex[atom->getLocalIndex()] = atom->getGlobalIndex() + 1;
723 }
724
725 //local index of cutoff group is trivial, it only depends on the order of travesing
726 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
727 localToGlobalCutoffGroupIndex.push_back(cg->getGlobalIndex() + 1);
728 }
729
730 }
731
732 //fill up mpiSimData struct
733 parallelData.nMolGlobal = getNGlobalMolecules();
734 parallelData.nMolLocal = getNMolecules();
735 parallelData.nAtomsGlobal = getNGlobalAtoms();
736 parallelData.nAtomsLocal = getNAtoms();
737 parallelData.nGroupsGlobal = getNGlobalCutoffGroups();
738 parallelData.nGroupsLocal = getNCutoffGroups();
739 parallelData.myNode = worldRank;
740 MPI_Comm_size(MPI_COMM_WORLD, &(parallelData.nProcessors));
741
742 //pass mpiSimData struct and index arrays to fortran
743 setFsimParallel(&parallelData, &(parallelData.nAtomsLocal),
744 &localToGlobalAtomIndex[0], &(parallelData.nGroupsLocal),
745 &localToGlobalCutoffGroupIndex[0], &isError);
746
747 if (isError) {
748 sprintf(painCave.errMsg,
749 "mpiRefresh errror: fortran didn't like something we gave it.\n");
750 painCave.isFatal = 1;
751 simError();
752 }
753
754 sprintf(checkPointMsg, " mpiRefresh successful.\n");
755 MPIcheckPoint();
756
757
758 }
759
760 #endif
761
762 double SimInfo::calcMaxCutoffRadius() {
763
764
765 std::set<AtomType*> atomTypes;
766 std::set<AtomType*>::iterator i;
767 std::vector<double> cutoffRadius;
768
769 //get the unique atom types
770 atomTypes = getUniqueAtomTypes();
771
772 //query the max cutoff radius among these atom types
773 for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
774 cutoffRadius.push_back(forceField_->getRcutFromAtomType(*i));
775 }
776
777 double maxCutoffRadius = *(std::max_element(cutoffRadius.begin(), cutoffRadius.end()));
778 #ifdef IS_MPI
779 //pick the max cutoff radius among the processors
780 #endif
781
782 return maxCutoffRadius;
783 }
784
785 void SimInfo::getCutoff(double& rcut, double& rsw) {
786
787 if (fInfo_.SIM_uses_Charges | fInfo_.SIM_uses_Dipoles | fInfo_.SIM_uses_RF) {
788
789 if (!simParams_->haveRcut()){
790 sprintf(painCave.errMsg,
791 "SimCreator Warning: No value was set for the cutoffRadius.\n"
792 "\tOOPSE will use a default value of 15.0 angstroms"
793 "\tfor the cutoffRadius.\n");
794 painCave.isFatal = 0;
795 simError();
796 rcut = 15.0;
797 } else{
798 rcut = simParams_->getRcut();
799 }
800
801 if (!simParams_->haveRsw()){
802 sprintf(painCave.errMsg,
803 "SimCreator Warning: No value was set for switchingRadius.\n"
804 "\tOOPSE will use a default value of\n"
805 "\t0.95 * cutoffRadius for the switchingRadius\n");
806 painCave.isFatal = 0;
807 simError();
808 rsw = 0.95 * rcut;
809 } else{
810 rsw = simParams_->getRsw();
811 }
812
813 } else {
814 // if charge, dipole or reaction field is not used and the cutofff radius is not specified in
815 //meta-data file, the maximum cutoff radius calculated from forcefiled will be used
816
817 if (simParams_->haveRcut()) {
818 rcut = simParams_->getRcut();
819 } else {
820 //set cutoff radius to the maximum cutoff radius based on atom types in the whole system
821 rcut = calcMaxCutoffRadius();
822 }
823
824 if (simParams_->haveRsw()) {
825 rsw = simParams_->getRsw();
826 } else {
827 rsw = rcut;
828 }
829
830 }
831 }
832
833 void SimInfo::setupCutoff() {
834 getCutoff(rcut_, rsw_);
835 double rnblist = rcut_ + 1; // skin of neighbor list
836
837 //Pass these cutoff radius etc. to fortran. This function should be called once and only once
838
839 int cp = TRADITIONAL_CUTOFF_POLICY;
840 if (simParams_->haveCutoffPolicy()) {
841 std::string myPolicy = simParams_->getCutoffPolicy();
842 if (myPolicy == "MIX") {
843 cp = MIX_CUTOFF_POLICY;
844 } else {
845 if (myPolicy == "MAX") {
846 cp = MAX_CUTOFF_POLICY;
847 } else {
848 if (myPolicy == "TRADITIONAL") {
849 cp = TRADITIONAL_CUTOFF_POLICY;
850 } else {
851 // throw error
852 sprintf( painCave.errMsg,
853 "SimInfo error: Unknown cutoffPolicy. (Input file specified %s .)\n\tcutoffPolicy must be one of: \"Mix\", \"Max\", or \"Traditional\".", myPolicy.c_str() );
854 painCave.isFatal = 1;
855 simError();
856 }
857 }
858 }
859 }
860 notifyFortranCutoffs(&rcut_, &rsw_, &rnblist, &cp);
861 // also send cutoff notification to electrostatics
862 setElectrostaticCutoffRadius(&rcut_);
863 }
864
865 void SimInfo::setupElectrostaticSummationMethod( int isError ) {
866
867 int errorOut;
868 int esm = NONE;
869 double alphaVal;
870 double dielectric;
871
872 errorOut = isError;
873 alphaVal = simParams_->getDampingAlpha();
874 dielectric = simParams_->getDielectric();
875
876 if (simParams_->haveElectrostaticSummationMethod()) {
877 std::string myMethod = simParams_->getElectrostaticSummationMethod();
878 if (myMethod == "NONE") {
879 esm = NONE;
880 } else {
881 if (myMethod == "UNDAMPED_WOLF") {
882 esm = UNDAMPED_WOLF;
883 } else {
884 if (myMethod == "DAMPED_WOLF") {
885 esm = DAMPED_WOLF;
886 if (!simParams_->haveDampingAlpha()) {
887 //throw error
888 sprintf( painCave.errMsg,
889 "SimInfo warning: dampingAlpha was not specified in the input file. A default value of %f (1/ang) will be used for the Damped Wolf Method.", alphaVal);
890 painCave.isFatal = 0;
891 simError();
892 }
893 } else {
894 if (myMethod == "REACTION_FIELD") {
895 esm = REACTION_FIELD;
896 } else {
897 // throw error
898 sprintf( painCave.errMsg,
899 "SimInfo error: Unknown electrostaticSummationMethod. (Input file specified %s .)\n\telectrostaticSummationMethod must be one of: \"none\", \"undamped_wolf\", \"damped_wolf\", or \"reaction_field\".", myMethod.c_str() );
900 painCave.isFatal = 1;
901 simError();
902 }
903 }
904 }
905 }
906 }
907 // let's pass some summation method variables to fortran
908 setElectrostaticSummationMethod( &esm );
909 setDampedWolfAlpha( &alphaVal );
910 setReactionFieldDielectric( &dielectric );
911 initFortranFF( &esm, &errorOut );
912 }
913
914 void SimInfo::addProperty(GenericData* genData) {
915 properties_.addProperty(genData);
916 }
917
918 void SimInfo::removeProperty(const std::string& propName) {
919 properties_.removeProperty(propName);
920 }
921
922 void SimInfo::clearProperties() {
923 properties_.clearProperties();
924 }
925
926 std::vector<std::string> SimInfo::getPropertyNames() {
927 return properties_.getPropertyNames();
928 }
929
930 std::vector<GenericData*> SimInfo::getProperties() {
931 return properties_.getProperties();
932 }
933
934 GenericData* SimInfo::getPropertyByName(const std::string& propName) {
935 return properties_.getPropertyByName(propName);
936 }
937
938 void SimInfo::setSnapshotManager(SnapshotManager* sman) {
939 if (sman_ == sman) {
940 return;
941 }
942 delete sman_;
943 sman_ = sman;
944
945 Molecule* mol;
946 RigidBody* rb;
947 Atom* atom;
948 SimInfo::MoleculeIterator mi;
949 Molecule::RigidBodyIterator rbIter;
950 Molecule::AtomIterator atomIter;;
951
952 for (mol = beginMolecule(mi); mol != NULL; mol = nextMolecule(mi)) {
953
954 for (atom = mol->beginAtom(atomIter); atom != NULL; atom = mol->nextAtom(atomIter)) {
955 atom->setSnapshotManager(sman_);
956 }
957
958 for (rb = mol->beginRigidBody(rbIter); rb != NULL; rb = mol->nextRigidBody(rbIter)) {
959 rb->setSnapshotManager(sman_);
960 }
961 }
962
963 }
964
965 Vector3d SimInfo::getComVel(){
966 SimInfo::MoleculeIterator i;
967 Molecule* mol;
968
969 Vector3d comVel(0.0);
970 double totalMass = 0.0;
971
972
973 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
974 double mass = mol->getMass();
975 totalMass += mass;
976 comVel += mass * mol->getComVel();
977 }
978
979 #ifdef IS_MPI
980 double tmpMass = totalMass;
981 Vector3d tmpComVel(comVel);
982 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
983 MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
984 #endif
985
986 comVel /= totalMass;
987
988 return comVel;
989 }
990
991 Vector3d SimInfo::getCom(){
992 SimInfo::MoleculeIterator i;
993 Molecule* mol;
994
995 Vector3d com(0.0);
996 double totalMass = 0.0;
997
998 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
999 double mass = mol->getMass();
1000 totalMass += mass;
1001 com += mass * mol->getCom();
1002 }
1003
1004 #ifdef IS_MPI
1005 double tmpMass = totalMass;
1006 Vector3d tmpCom(com);
1007 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1008 MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1009 #endif
1010
1011 com /= totalMass;
1012
1013 return com;
1014
1015 }
1016
1017 std::ostream& operator <<(std::ostream& o, SimInfo& info) {
1018
1019 return o;
1020 }
1021
1022
1023 /*
1024 Returns center of mass and center of mass velocity in one function call.
1025 */
1026
1027 void SimInfo::getComAll(Vector3d &com, Vector3d &comVel){
1028 SimInfo::MoleculeIterator i;
1029 Molecule* mol;
1030
1031
1032 double totalMass = 0.0;
1033
1034
1035 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1036 double mass = mol->getMass();
1037 totalMass += mass;
1038 com += mass * mol->getCom();
1039 comVel += mass * mol->getComVel();
1040 }
1041
1042 #ifdef IS_MPI
1043 double tmpMass = totalMass;
1044 Vector3d tmpCom(com);
1045 Vector3d tmpComVel(comVel);
1046 MPI_Allreduce(&tmpMass,&totalMass,1,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1047 MPI_Allreduce(tmpCom.getArrayPointer(), com.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1048 MPI_Allreduce(tmpComVel.getArrayPointer(), comVel.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1049 #endif
1050
1051 com /= totalMass;
1052 comVel /= totalMass;
1053 }
1054
1055 /*
1056 Return intertia tensor for entire system and angular momentum Vector.
1057
1058
1059 [ Ixx -Ixy -Ixz ]
1060 J =| -Iyx Iyy -Iyz |
1061 [ -Izx -Iyz Izz ]
1062 */
1063
1064 void SimInfo::getInertiaTensor(Mat3x3d &inertiaTensor, Vector3d &angularMomentum){
1065
1066
1067 double xx = 0.0;
1068 double yy = 0.0;
1069 double zz = 0.0;
1070 double xy = 0.0;
1071 double xz = 0.0;
1072 double yz = 0.0;
1073 Vector3d com(0.0);
1074 Vector3d comVel(0.0);
1075
1076 getComAll(com, comVel);
1077
1078 SimInfo::MoleculeIterator i;
1079 Molecule* mol;
1080
1081 Vector3d thisq(0.0);
1082 Vector3d thisv(0.0);
1083
1084 double thisMass = 0.0;
1085
1086
1087
1088
1089 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1090
1091 thisq = mol->getCom()-com;
1092 thisv = mol->getComVel()-comVel;
1093 thisMass = mol->getMass();
1094 // Compute moment of intertia coefficients.
1095 xx += thisq[0]*thisq[0]*thisMass;
1096 yy += thisq[1]*thisq[1]*thisMass;
1097 zz += thisq[2]*thisq[2]*thisMass;
1098
1099 // compute products of intertia
1100 xy += thisq[0]*thisq[1]*thisMass;
1101 xz += thisq[0]*thisq[2]*thisMass;
1102 yz += thisq[1]*thisq[2]*thisMass;
1103
1104 angularMomentum += cross( thisq, thisv ) * thisMass;
1105
1106 }
1107
1108
1109 inertiaTensor(0,0) = yy + zz;
1110 inertiaTensor(0,1) = -xy;
1111 inertiaTensor(0,2) = -xz;
1112 inertiaTensor(1,0) = -xy;
1113 inertiaTensor(1,1) = xx + zz;
1114 inertiaTensor(1,2) = -yz;
1115 inertiaTensor(2,0) = -xz;
1116 inertiaTensor(2,1) = -yz;
1117 inertiaTensor(2,2) = xx + yy;
1118
1119 #ifdef IS_MPI
1120 Mat3x3d tmpI(inertiaTensor);
1121 Vector3d tmpAngMom;
1122 MPI_Allreduce(tmpI.getArrayPointer(), inertiaTensor.getArrayPointer(),9,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1123 MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1124 #endif
1125
1126 return;
1127 }
1128
1129 //Returns the angular momentum of the system
1130 Vector3d SimInfo::getAngularMomentum(){
1131
1132 Vector3d com(0.0);
1133 Vector3d comVel(0.0);
1134 Vector3d angularMomentum(0.0);
1135
1136 getComAll(com,comVel);
1137
1138 SimInfo::MoleculeIterator i;
1139 Molecule* mol;
1140
1141 Vector3d thisr(0.0);
1142 Vector3d thisp(0.0);
1143
1144 double thisMass;
1145
1146 for (mol = beginMolecule(i); mol != NULL; mol = nextMolecule(i)) {
1147 thisMass = mol->getMass();
1148 thisr = mol->getCom()-com;
1149 thisp = (mol->getComVel()-comVel)*thisMass;
1150
1151 angularMomentum += cross( thisr, thisp );
1152
1153 }
1154
1155 #ifdef IS_MPI
1156 Vector3d tmpAngMom;
1157 MPI_Allreduce(tmpAngMom.getArrayPointer(), angularMomentum.getArrayPointer(),3,MPI_DOUBLE,MPI_SUM, MPI_COMM_WORLD);
1158 #endif
1159
1160 return angularMomentum;
1161 }
1162
1163
1164 }//end namespace oopse
1165