ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/brains/SimInfo.cpp
Revision: 2285
Committed: Wed Sep 7 20:46:46 2005 UTC (18 years, 10 months ago) by gezelter
File size: 33570 byte(s)
Log Message:
adding c-side interface to change cutoff Policy

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