ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/ForceManager.cpp
Revision: 1545
Committed: Fri Apr 8 21:25:19 2011 UTC (14 years ago) by gezelter
File size: 18356 byte(s)
Log Message:
still busted, but much progress

File Contents

# User Rev Content
1 gezelter 507 /*
2 gezelter 246 * 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 gezelter 1390 * 1. Redistributions of source code must retain the above copyright
10 gezelter 246 * notice, this list of conditions and the following disclaimer.
11     *
12 gezelter 1390 * 2. Redistributions in binary form must reproduce the above copyright
13 gezelter 246 * notice, this list of conditions and the following disclaimer in the
14     * documentation and/or other materials provided with the
15     * distribution.
16     *
17     * This software is provided "AS IS," without a warranty of any
18     * kind. All express or implied conditions, representations and
19     * warranties, including any implied warranty of merchantability,
20     * fitness for a particular purpose or non-infringement, are hereby
21     * excluded. The University of Notre Dame and its licensors shall not
22     * be liable for any damages suffered by licensee as a result of
23     * using, modifying or distributing the software or its
24     * derivatives. In no event will the University of Notre Dame or its
25     * licensors be liable for any lost revenue, profit or data, or for
26     * direct, indirect, special, consequential, incidental or punitive
27     * damages, however caused and regardless of the theory of liability,
28     * arising out of the use of or inability to use software, even if the
29     * University of Notre Dame has been advised of the possibility of
30     * such damages.
31 gezelter 1390 *
32     * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33     * research, please cite the appropriate papers when you publish your
34     * work. Good starting points are:
35     *
36     * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37     * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38     * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).
39     * [4] Vardeman & Gezelter, in progress (2009).
40 gezelter 246 */
41    
42 gezelter 507 /**
43     * @file ForceManager.cpp
44     * @author tlin
45     * @date 11/09/2004
46     * @time 10:39am
47     * @version 1.0
48     */
49 gezelter 246
50     #include "brains/ForceManager.hpp"
51     #include "primitives/Molecule.hpp"
52     #include "UseTheForce/doForces_interface.h"
53 gezelter 1390 #define __OPENMD_C
54 chuckv 664 #include "UseTheForce/DarkSide/fInteractionMap.h"
55 gezelter 246 #include "utils/simError.h"
56 xsun 1215 #include "primitives/Bond.hpp"
57 tim 749 #include "primitives/Bend.hpp"
58 cli2 1275 #include "primitives/Torsion.hpp"
59     #include "primitives/Inversion.hpp"
60 gezelter 1544 #include "parallel/ForceDecomposition.hpp"
61     //#include "parallel/SerialDecomposition.hpp"
62 gezelter 1467
63 gezelter 1545 using namespace std;
64 gezelter 1390 namespace OpenMD {
65 gezelter 1469
66 gezelter 1545 ForceManager::ForceManager(SimInfo * info) : info_(info) {
67    
68 gezelter 1544 #ifdef IS_MPI
69     decomp_ = new ForceDecomposition(info_);
70     #else
71 gezelter 1545 // decomp_ = new SerialDecomposition(info);
72 gezelter 1544 #endif
73 gezelter 1469 }
74 gezelter 1545
75 gezelter 1464 void ForceManager::calcForces() {
76 gezelter 1126
77 gezelter 246 if (!info_->isFortranInitialized()) {
78 gezelter 507 info_->update();
79 gezelter 1535 nbiMan_->setSimInfo(info_);
80 gezelter 1544 nbiMan_->initialize();
81 gezelter 1545 swfun_ = nbiMan_->getSwitchingFunction();
82 gezelter 1544 decomp_->distributeInitialData();
83 gezelter 1535 info_->setupFortran();
84 gezelter 246 }
85 gezelter 1126
86 gezelter 1544 preCalculation();
87 gezelter 246 calcShortRangeInteraction();
88 gezelter 1464 calcLongRangeInteraction();
89     postCalculation();
90 tim 749
91 gezelter 507 }
92 gezelter 1126
93 gezelter 507 void ForceManager::preCalculation() {
94 gezelter 246 SimInfo::MoleculeIterator mi;
95     Molecule* mol;
96     Molecule::AtomIterator ai;
97     Atom* atom;
98     Molecule::RigidBodyIterator rbIter;
99     RigidBody* rb;
100 gezelter 1540 Molecule::CutoffGroupIterator ci;
101     CutoffGroup* cg;
102 gezelter 246
103     // forces are zeroed here, before any are accumulated.
104 chuckv 1245
105 gezelter 1126 for (mol = info_->beginMolecule(mi); mol != NULL;
106     mol = info_->nextMolecule(mi)) {
107 gezelter 507 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
108     atom->zeroForcesAndTorques();
109     }
110 chuckv 1245
111 gezelter 507 //change the positions of atoms which belong to the rigidbodies
112 gezelter 1126 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
113     rb = mol->nextRigidBody(rbIter)) {
114 gezelter 507 rb->zeroForcesAndTorques();
115     }
116 gezelter 1540
117     if(info_->getNGlobalCutoffGroups() != info_->getNGlobalAtoms()){
118     for(cg = mol->beginCutoffGroup(ci); cg != NULL;
119     cg = mol->nextCutoffGroup(ci)) {
120     //calculate the center of mass of cutoff group
121     cg->updateCOM();
122     }
123     }
124 gezelter 246 }
125 gezelter 1540
126 gezelter 1126 // Zero out the stress tensor
127     tau *= 0.0;
128    
129 gezelter 507 }
130 gezelter 1126
131 gezelter 507 void ForceManager::calcShortRangeInteraction() {
132 gezelter 246 Molecule* mol;
133     RigidBody* rb;
134     Bond* bond;
135     Bend* bend;
136     Torsion* torsion;
137 cli2 1275 Inversion* inversion;
138 gezelter 246 SimInfo::MoleculeIterator mi;
139     Molecule::RigidBodyIterator rbIter;
140     Molecule::BondIterator bondIter;;
141     Molecule::BendIterator bendIter;
142     Molecule::TorsionIterator torsionIter;
143 cli2 1275 Molecule::InversionIterator inversionIter;
144 tim 963 RealType bondPotential = 0.0;
145     RealType bendPotential = 0.0;
146     RealType torsionPotential = 0.0;
147 cli2 1275 RealType inversionPotential = 0.0;
148 gezelter 246
149     //calculate short range interactions
150 gezelter 1126 for (mol = info_->beginMolecule(mi); mol != NULL;
151     mol = info_->nextMolecule(mi)) {
152 gezelter 246
153 gezelter 507 //change the positions of atoms which belong to the rigidbodies
154 gezelter 1126 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
155     rb = mol->nextRigidBody(rbIter)) {
156     rb->updateAtoms();
157 gezelter 507 }
158 gezelter 246
159 gezelter 1126 for (bond = mol->beginBond(bondIter); bond != NULL;
160     bond = mol->nextBond(bondIter)) {
161 tim 749 bond->calcForce();
162     bondPotential += bond->getPotential();
163 gezelter 507 }
164 gezelter 246
165 gezelter 1126 for (bend = mol->beginBend(bendIter); bend != NULL;
166     bend = mol->nextBend(bendIter)) {
167    
168     RealType angle;
169     bend->calcForce(angle);
170     RealType currBendPot = bend->getPotential();
171 gezelter 1448
172 gezelter 1126 bendPotential += bend->getPotential();
173 gezelter 1545 map<Bend*, BendDataSet>::iterator i = bendDataSets.find(bend);
174 gezelter 1126 if (i == bendDataSets.end()) {
175     BendDataSet dataSet;
176     dataSet.prev.angle = dataSet.curr.angle = angle;
177     dataSet.prev.potential = dataSet.curr.potential = currBendPot;
178     dataSet.deltaV = 0.0;
179 gezelter 1545 bendDataSets.insert(map<Bend*, BendDataSet>::value_type(bend, dataSet));
180 gezelter 1126 }else {
181     i->second.prev.angle = i->second.curr.angle;
182     i->second.prev.potential = i->second.curr.potential;
183     i->second.curr.angle = angle;
184     i->second.curr.potential = currBendPot;
185     i->second.deltaV = fabs(i->second.curr.potential -
186     i->second.prev.potential);
187     }
188 gezelter 507 }
189 gezelter 1126
190     for (torsion = mol->beginTorsion(torsionIter); torsion != NULL;
191     torsion = mol->nextTorsion(torsionIter)) {
192 tim 963 RealType angle;
193 gezelter 1126 torsion->calcForce(angle);
194 tim 963 RealType currTorsionPot = torsion->getPotential();
195 gezelter 1126 torsionPotential += torsion->getPotential();
196 gezelter 1545 map<Torsion*, TorsionDataSet>::iterator i = torsionDataSets.find(torsion);
197 gezelter 1126 if (i == torsionDataSets.end()) {
198     TorsionDataSet dataSet;
199     dataSet.prev.angle = dataSet.curr.angle = angle;
200     dataSet.prev.potential = dataSet.curr.potential = currTorsionPot;
201     dataSet.deltaV = 0.0;
202 gezelter 1545 torsionDataSets.insert(map<Torsion*, TorsionDataSet>::value_type(torsion, dataSet));
203 gezelter 1126 }else {
204     i->second.prev.angle = i->second.curr.angle;
205     i->second.prev.potential = i->second.curr.potential;
206     i->second.curr.angle = angle;
207     i->second.curr.potential = currTorsionPot;
208     i->second.deltaV = fabs(i->second.curr.potential -
209     i->second.prev.potential);
210     }
211     }
212 gezelter 1545
213 cli2 1275 for (inversion = mol->beginInversion(inversionIter);
214     inversion != NULL;
215     inversion = mol->nextInversion(inversionIter)) {
216     RealType angle;
217     inversion->calcForce(angle);
218     RealType currInversionPot = inversion->getPotential();
219     inversionPotential += inversion->getPotential();
220 gezelter 1545 map<Inversion*, InversionDataSet>::iterator i = inversionDataSets.find(inversion);
221 cli2 1275 if (i == inversionDataSets.end()) {
222     InversionDataSet dataSet;
223     dataSet.prev.angle = dataSet.curr.angle = angle;
224     dataSet.prev.potential = dataSet.curr.potential = currInversionPot;
225     dataSet.deltaV = 0.0;
226 gezelter 1545 inversionDataSets.insert(map<Inversion*, InversionDataSet>::value_type(inversion, dataSet));
227 cli2 1275 }else {
228     i->second.prev.angle = i->second.curr.angle;
229     i->second.prev.potential = i->second.curr.potential;
230     i->second.curr.angle = angle;
231     i->second.curr.potential = currInversionPot;
232     i->second.deltaV = fabs(i->second.curr.potential -
233     i->second.prev.potential);
234     }
235     }
236 gezelter 246 }
237    
238 gezelter 1126 RealType shortRangePotential = bondPotential + bendPotential +
239 cli2 1275 torsionPotential + inversionPotential;
240 gezelter 246 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
241     curSnapshot->statData[Stats::SHORT_RANGE_POTENTIAL] = shortRangePotential;
242 tim 665 curSnapshot->statData[Stats::BOND_POTENTIAL] = bondPotential;
243     curSnapshot->statData[Stats::BEND_POTENTIAL] = bendPotential;
244     curSnapshot->statData[Stats::DIHEDRAL_POTENTIAL] = torsionPotential;
245 gezelter 1545 curSnapshot->statData[Stats::INVERSION_POTENTIAL] = inversionPotential;
246 gezelter 507 }
247 gezelter 1126
248 gezelter 1464 void ForceManager::calcLongRangeInteraction() {
249 gezelter 246
250 gezelter 1545 // some of this initial stuff will go away:
251     Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
252     DataStorage* config = &(curSnapshot->atomData);
253     DataStorage* cgConfig = &(curSnapshot->cgData);
254     RealType* frc = config->getArrayPointer(DataStorage::dslForce);
255     RealType* pos = config->getArrayPointer(DataStorage::dslPosition);
256     RealType* trq = config->getArrayPointer(DataStorage::dslTorque);
257     RealType* A = config->getArrayPointer(DataStorage::dslAmat);
258     RealType* electroFrame = config->getArrayPointer(DataStorage::dslElectroFrame);
259     RealType* particlePot = config->getArrayPointer(DataStorage::dslParticlePot);
260     RealType* rc;
261    
262 gezelter 1540 if(info_->getNGlobalCutoffGroups() != info_->getNGlobalAtoms()){
263     rc = cgConfig->getArrayPointer(DataStorage::dslPosition);
264 gezelter 246 } else {
265 gezelter 1126 // center of mass of the group is the same as position of the atom
266     // if cutoff group does not exist
267 gezelter 507 rc = pos;
268 gezelter 246 }
269 gezelter 1126
270 gezelter 246 //initialize data before passing to fortran
271 tim 963 RealType longRangePotential[LR_POT_TYPES];
272     RealType lrPot = 0.0;
273 gezelter 246 int isError = 0;
274    
275 chuckv 664 for (int i=0; i<LR_POT_TYPES;i++){
276     longRangePotential[i]=0.0; //Initialize array
277     }
278 gezelter 1545
279     // new stuff starts here:
280    
281 gezelter 1544 decomp_->distributeData();
282 gezelter 1545
283     int cg1, cg2;
284     Vector3d d_grp;
285     RealType rgrpsq, rgrp;
286     RealType vij;
287     Vector3d fij, fg;
288     pair<int, int> gtypes;
289     RealType rCutSq;
290     bool in_switching_region;
291     RealType sw, dswdr, swderiv;
292     vector<int> atomListI;
293     vector<int> atomListJ;
294     InteractionData idat;
295 gezelter 1544
296 gezelter 1545 int loopStart, loopEnd;
297 gezelter 1544
298 gezelter 1545 loopEnd = PAIR_LOOP;
299     if (info_->requiresPrepair_) {
300     loopStart = PREPAIR_LOOP;
301     } else {
302     loopStart = PAIR_LOOP;
303     }
304    
305     for (int iLoop = loopStart; iLoop < loopEnd; iLoop++) {
306    
307     if (iLoop == loopStart) {
308     bool update_nlist = decomp_->checkNeighborList();
309     if (update_nlist)
310     neighborList = decomp_->buildNeighborList();
311 gezelter 1544 }
312 gezelter 1545
313     for (vector<pair<int, int> >::iterator it = neighborList.begin();
314     it != neighborList.end(); ++it) {
315    
316     cg1 = (*it).first;
317     cg2 = (*it).second;
318    
319     gtypes = decomp_->getGroupTypes(cg1, cg2);
320     d_grp = decomp_->getIntergroupVector(cg1, cg2);
321     curSnapshot->wrapVector(d_grp);
322     rgrpsq = d_grp.lengthSquare();
323     rCutSq = groupCutoffMap(gtypes).first;
324    
325     if (rgrpsq < rCutSq) {
326     idat.rcut = groupCutoffMap(gtypes).second;
327     if (iLoop == PAIR_LOOP) {
328     vij = 0.0;
329     fij = V3Zero;
330     }
331    
332     in_switching_region = swfun_->getSwitch(rgrpsq, idat.sw, idat.dswdr, rgrp);
333    
334     atomListI = decomp_->getAtomsInGroupI(cg1);
335     atomListJ = decomp_->getAtomsInGroupJ(cg2);
336    
337     for (vector<int>::iterator ia = atomListI.begin();
338     ia != atomListI.end(); ++ia) {
339     atom1 = (*ia);
340    
341     for (vector<int>::iterator jb = atomListJ.begin();
342     jb != atomListJ.end(); ++jb) {
343     atom2 = (*jb);
344    
345     if (!decomp_->skipAtomPair(atom1, atom2)) {
346    
347     if (atomListI.size() == 1 && atomListJ.size() == 1) {
348     idat.d = d_grp;
349     idat.r2 = rgrpsq;
350     } else {
351     idat.d = decomp_->getInteratomicVector(atom1, atom2);
352     curSnapshot->wrapVector(idat.d);
353     idat.r2 = idat.d.lengthSquare();
354     }
355    
356     idat.r = sqrt(idat.r2);
357     decomp_->fillInteractionData(atom1, atom2, idat);
358    
359     if (iLoop == PREPAIR_LOOP) {
360     interactionMan_->doPrePair(idat);
361     } else {
362     interactionMan_->doPair(idat);
363     vij += idat.vpair;
364     fij += idat.f1;
365     tau -= outProduct(idat.d, idat.f);
366     }
367     }
368     }
369     }
370    
371     if (iLoop == PAIR_LOOP) {
372     if (in_switching_region) {
373     swderiv = vij * dswdr / rgrp;
374     fg = swderiv * d_grp;
375    
376     fij += fg;
377    
378     if (atomListI.size() == 1 && atomListJ.size() == 1) {
379     tau -= outProduct(idat.d, fg);
380     }
381    
382     for (vector<int>::iterator ia = atomListI.begin();
383     ia != atomListI.end(); ++ia) {
384     atom1 = (*ia);
385     mf = decomp_->getMfactI(atom1);
386     // fg is the force on atom ia due to cutoff group's
387     // presence in switching region
388     fg = swderiv * d_grp * mf;
389     decomp_->addForceToAtomI(atom1, fg);
390    
391     if (atomListI.size() > 1) {
392     if (info_->usesAtomicVirial_) {
393     // find the distance between the atom
394     // and the center of the cutoff group:
395     dag = decomp_->getAtomToGroupVectorI(atom1, cg1);
396     tau -= outProduct(dag, fg);
397     }
398     }
399     }
400     for (vector<int>::iterator jb = atomListJ.begin();
401     jb != atomListJ.end(); ++jb) {
402     atom2 = (*jb);
403     mf = decomp_->getMfactJ(atom2);
404     // fg is the force on atom jb due to cutoff group's
405     // presence in switching region
406     fg = -swderiv * d_grp * mf;
407     decomp_->addForceToAtomJ(atom2, fg);
408    
409     if (atomListJ.size() > 1) {
410     if (info_->usesAtomicVirial_) {
411     // find the distance between the atom
412     // and the center of the cutoff group:
413     dag = decomp_->getAtomToGroupVectorJ(atom2, cg2);
414     tau -= outProduct(dag, fg);
415     }
416     }
417     }
418     }
419     //if (!SIM_uses_AtomicVirial) {
420     // tau -= outProduct(d_grp, fij);
421     //}
422     }
423     }
424     }
425    
426     if (iLoop == PREPAIR_LOOP) {
427     if (info_->requiresPrepair_) {
428     decomp_->collectIntermediateData();
429     atomList = decomp_->getAtomList();
430     for (vector<int>::iterator ia = atomList.begin();
431     ia != atomList.end(); ++ia) {
432     atom1 = (*ia);
433     decomp_->populateSelfData(atom1, SelfData sdat);
434     interactionMan_->doPreForce(sdat);
435     }
436     decomp_->distributeIntermediateData();
437     }
438     }
439    
440 gezelter 1544 }
441 gezelter 1545
442 gezelter 1544 decomp_->collectData();
443 gezelter 1545
444     if (info_->requiresSkipCorrection_ || info_->requiresSelfCorrection_) {
445     atomList = decomp_->getAtomList();
446     for (vector<int>::iterator ia = atomList.begin();
447     ia != atomList.end(); ++ia) {
448     atom1 = (*ia);
449 gezelter 1544
450 gezelter 1545 if (info_->requiresSkipCorrection_) {
451     vector<int> skipList = decomp_->getSkipsForAtom(atom1);
452     for (vector<int>::iterator jb = skipList.begin();
453     jb != skipList.end(); ++jb) {
454     atom2 = (*jb);
455     decomp_->populateSkipData(atom1, atom2, InteractionData idat);
456     interactionMan_->doSkipCorrection(idat);
457     }
458     }
459    
460     if (info_->requiresSelfCorrection_) {
461     decomp_->populateSelfData(atom1, SelfData sdat);
462     interactionMan_->doSelfCorrection(sdat);
463     }
464    
465    
466 gezelter 246 }
467 gezelter 1545
468 chuckv 664 for (int i=0; i<LR_POT_TYPES;i++){
469     lrPot += longRangePotential[i]; //Quick hack
470     }
471 gezelter 1503
472 gezelter 246 //store the tau and long range potential
473 chuckv 664 curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL] = lrPot;
474 chrisfen 691 curSnapshot->statData[Stats::VANDERWAALS_POTENTIAL] = longRangePotential[VDW_POT];
475 tim 681 curSnapshot->statData[Stats::ELECTROSTATIC_POTENTIAL] = longRangePotential[ELECTROSTATIC_POT];
476 gezelter 507 }
477 gezelter 246
478 gezelter 1126
479 gezelter 1464 void ForceManager::postCalculation() {
480 gezelter 246 SimInfo::MoleculeIterator mi;
481     Molecule* mol;
482     Molecule::RigidBodyIterator rbIter;
483     RigidBody* rb;
484 gezelter 1126 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
485 gezelter 246
486     // collect the atomic forces onto rigid bodies
487 gezelter 1126
488     for (mol = info_->beginMolecule(mi); mol != NULL;
489     mol = info_->nextMolecule(mi)) {
490     for (rb = mol->beginRigidBody(rbIter); rb != NULL;
491     rb = mol->nextRigidBody(rbIter)) {
492 gezelter 1464 Mat3x3d rbTau = rb->calcForcesAndTorquesAndVirial();
493     tau += rbTau;
494 gezelter 507 }
495 gezelter 1126 }
496 gezelter 1464
497 gezelter 1126 #ifdef IS_MPI
498 gezelter 1464 Mat3x3d tmpTau(tau);
499     MPI_Allreduce(tmpTau.getArrayPointer(), tau.getArrayPointer(),
500     9, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
501 gezelter 1126 #endif
502 gezelter 1464 curSnapshot->statData.setTau(tau);
503 gezelter 507 }
504 gezelter 246
505 gezelter 1390 } //end namespace OpenMD

Properties

Name Value
svn:executable *
svn:keywords Author Id Revision Date