ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/ForceManager.cpp
Revision: 1581
Committed: Mon Jun 13 22:13:12 2011 UTC (13 years, 10 months ago) by gezelter
File size: 26997 byte(s)
Log Message:
bug fixes

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 gezelter 1576
51 gezelter 246 #include "brains/ForceManager.hpp"
52     #include "primitives/Molecule.hpp"
53 gezelter 1390 #define __OPENMD_C
54 gezelter 246 #include "utils/simError.h"
55 xsun 1215 #include "primitives/Bond.hpp"
56 tim 749 #include "primitives/Bend.hpp"
57 cli2 1275 #include "primitives/Torsion.hpp"
58     #include "primitives/Inversion.hpp"
59 gezelter 1551 #include "nonbonded/NonBondedInteraction.hpp"
60 gezelter 1549 #include "parallel/ForceMatrixDecomposition.hpp"
61 gezelter 1467
62 gezelter 1545 using namespace std;
63 gezelter 1390 namespace OpenMD {
64 gezelter 1469
65 gezelter 1545 ForceManager::ForceManager(SimInfo * info) : info_(info) {
66 gezelter 1576 forceField_ = info_->getForceField();
67 gezelter 1577 interactionMan_ = new InteractionManager();
68 gezelter 1579 fDecomp_ = new ForceMatrixDecomposition(info_, interactionMan_);
69 gezelter 1469 }
70 gezelter 1576
71     /**
72     * setupCutoffs
73     *
74     * Sets the values of cutoffRadius, cutoffMethod, and cutoffPolicy
75     *
76     * cutoffRadius : realType
77     * If the cutoffRadius was explicitly set, use that value.
78     * If the cutoffRadius was not explicitly set:
79     * Are there electrostatic atoms? Use 12.0 Angstroms.
80     * No electrostatic atoms? Poll the atom types present in the
81     * simulation for suggested cutoff values (e.g. 2.5 * sigma).
82     * Use the maximum suggested value that was found.
83     *
84     * cutoffMethod : (one of HARD, SWITCHED, SHIFTED_FORCE, SHIFTED_POTENTIAL)
85     * If cutoffMethod was explicitly set, use that choice.
86     * If cutoffMethod was not explicitly set, use SHIFTED_FORCE
87     *
88     * cutoffPolicy : (one of MIX, MAX, TRADITIONAL)
89     * If cutoffPolicy was explicitly set, use that choice.
90     * If cutoffPolicy was not explicitly set, use TRADITIONAL
91     */
92     void ForceManager::setupCutoffs() {
93    
94     Globals* simParams_ = info_->getSimParams();
95     ForceFieldOptions& forceFieldOptions_ = forceField_->getForceFieldOptions();
96    
97     if (simParams_->haveCutoffRadius()) {
98     rCut_ = simParams_->getCutoffRadius();
99     } else {
100     if (info_->usesElectrostaticAtoms()) {
101     sprintf(painCave.errMsg,
102     "ForceManager::setupCutoffs: No value was set for the cutoffRadius.\n"
103     "\tOpenMD will use a default value of 12.0 angstroms"
104     "\tfor the cutoffRadius.\n");
105     painCave.isFatal = 0;
106     painCave.severity = OPENMD_INFO;
107     simError();
108     rCut_ = 12.0;
109     } else {
110     RealType thisCut;
111     set<AtomType*>::iterator i;
112     set<AtomType*> atomTypes;
113     atomTypes = info_->getSimulatedAtomTypes();
114     for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
115     thisCut = interactionMan_->getSuggestedCutoffRadius((*i));
116     rCut_ = max(thisCut, rCut_);
117     }
118     sprintf(painCave.errMsg,
119     "ForceManager::setupCutoffs: No value was set for the cutoffRadius.\n"
120     "\tOpenMD will use %lf angstroms.\n",
121     rCut_);
122     painCave.isFatal = 0;
123     painCave.severity = OPENMD_INFO;
124     simError();
125 gezelter 1579 }
126     fDecomp_->setUserCutoff(rCut_);
127 gezelter 1576 }
128    
129     map<string, CutoffMethod> stringToCutoffMethod;
130     stringToCutoffMethod["HARD"] = HARD;
131     stringToCutoffMethod["SWITCHED"] = SWITCHED;
132     stringToCutoffMethod["SHIFTED_POTENTIAL"] = SHIFTED_POTENTIAL;
133     stringToCutoffMethod["SHIFTED_FORCE"] = SHIFTED_FORCE;
134 gezelter 1545
135 gezelter 1576 if (simParams_->haveCutoffMethod()) {
136     string cutMeth = toUpperCopy(simParams_->getCutoffMethod());
137     map<string, CutoffMethod>::iterator i;
138     i = stringToCutoffMethod.find(cutMeth);
139     if (i == stringToCutoffMethod.end()) {
140     sprintf(painCave.errMsg,
141     "ForceManager::setupCutoffs: Could not find chosen cutoffMethod %s\n"
142     "\tShould be one of: "
143     "HARD, SWITCHED, SHIFTED_POTENTIAL, or SHIFTED_FORCE\n",
144     cutMeth.c_str());
145     painCave.isFatal = 1;
146     painCave.severity = OPENMD_ERROR;
147     simError();
148     } else {
149     cutoffMethod_ = i->second;
150     }
151     } else {
152     sprintf(painCave.errMsg,
153     "ForceManager::setupCutoffs: No value was set for the cutoffMethod.\n"
154     "\tOpenMD will use SHIFTED_FORCE.\n");
155     painCave.isFatal = 0;
156     painCave.severity = OPENMD_INFO;
157     simError();
158     cutoffMethod_ = SHIFTED_FORCE;
159     }
160    
161     map<string, CutoffPolicy> stringToCutoffPolicy;
162     stringToCutoffPolicy["MIX"] = MIX;
163     stringToCutoffPolicy["MAX"] = MAX;
164     stringToCutoffPolicy["TRADITIONAL"] = TRADITIONAL;
165    
166     std::string cutPolicy;
167     if (forceFieldOptions_.haveCutoffPolicy()){
168     cutPolicy = forceFieldOptions_.getCutoffPolicy();
169     }else if (simParams_->haveCutoffPolicy()) {
170     cutPolicy = simParams_->getCutoffPolicy();
171     }
172    
173     if (!cutPolicy.empty()){
174     toUpper(cutPolicy);
175     map<string, CutoffPolicy>::iterator i;
176     i = stringToCutoffPolicy.find(cutPolicy);
177    
178     if (i == stringToCutoffPolicy.end()) {
179     sprintf(painCave.errMsg,
180     "ForceManager::setupCutoffs: Could not find chosen cutoffPolicy %s\n"
181     "\tShould be one of: "
182     "MIX, MAX, or TRADITIONAL\n",
183     cutPolicy.c_str());
184     painCave.isFatal = 1;
185     painCave.severity = OPENMD_ERROR;
186     simError();
187     } else {
188     cutoffPolicy_ = i->second;
189     }
190     } else {
191     sprintf(painCave.errMsg,
192     "ForceManager::setupCutoffs: No value was set for the cutoffPolicy.\n"
193     "\tOpenMD will use TRADITIONAL.\n");
194     painCave.isFatal = 0;
195     painCave.severity = OPENMD_INFO;
196     simError();
197     cutoffPolicy_ = TRADITIONAL;
198     }
199 gezelter 1579 fDecomp_->setCutoffPolicy(cutoffPolicy_);
200 gezelter 1576 }
201    
202     /**
203     * setupSwitching
204     *
205     * Sets the values of switchingRadius and
206     * If the switchingRadius was explicitly set, use that value (but check it)
207     * If the switchingRadius was not explicitly set: use 0.85 * cutoffRadius_
208     */
209     void ForceManager::setupSwitching() {
210     Globals* simParams_ = info_->getSimParams();
211 gezelter 1577
212     // create the switching function object:
213     switcher_ = new SwitchingFunction();
214 gezelter 1126
215 gezelter 1576 if (simParams_->haveSwitchingRadius()) {
216     rSwitch_ = simParams_->getSwitchingRadius();
217     if (rSwitch_ > rCut_) {
218     sprintf(painCave.errMsg,
219 gezelter 1577 "ForceManager::setupSwitching: switchingRadius (%f) is larger "
220     "than the cutoffRadius(%f)\n", rSwitch_, rCut_);
221 gezelter 1576 painCave.isFatal = 1;
222     painCave.severity = OPENMD_ERROR;
223     simError();
224     }
225     } else {
226     rSwitch_ = 0.85 * rCut_;
227     sprintf(painCave.errMsg,
228     "ForceManager::setupSwitching: No value was set for the switchingRadius.\n"
229     "\tOpenMD will use a default value of 85 percent of the cutoffRadius.\n"
230     "\tswitchingRadius = %f. for this simulation\n", rSwitch_);
231     painCave.isFatal = 0;
232     painCave.severity = OPENMD_WARNING;
233     simError();
234     }
235    
236 gezelter 1577 // Default to cubic switching function.
237     sft_ = cubic;
238 gezelter 1576 if (simParams_->haveSwitchingFunctionType()) {
239     string funcType = simParams_->getSwitchingFunctionType();
240     toUpper(funcType);
241     if (funcType == "CUBIC") {
242     sft_ = cubic;
243     } else {
244     if (funcType == "FIFTH_ORDER_POLYNOMIAL") {
245     sft_ = fifth_order_poly;
246     } else {
247     // throw error
248     sprintf( painCave.errMsg,
249     "ForceManager::setupSwitching : Unknown switchingFunctionType. (Input file specified %s .)\n"
250     "\tswitchingFunctionType must be one of: "
251     "\"cubic\" or \"fifth_order_polynomial\".",
252     funcType.c_str() );
253     painCave.isFatal = 1;
254     painCave.severity = OPENMD_ERROR;
255     simError();
256     }
257     }
258     }
259     switcher_->setSwitchType(sft_);
260     switcher_->setSwitch(rSwitch_, rCut_);
261     }
262    
263     void ForceManager::initialize() {
264    
265 gezelter 1569 if (!info_->isTopologyDone()) {
266 gezelter 507 info_->update();
267 gezelter 1546 interactionMan_->setSimInfo(info_);
268     interactionMan_->initialize();
269 gezelter 1576
270     // We want to delay the cutoffs until after the interaction
271     // manager has set up the atom-atom interactions so that we can
272     // query them for suggested cutoff values
273    
274     setupCutoffs();
275     setupSwitching();
276    
277     info_->prepareTopology();
278 gezelter 246 }
279 gezelter 1576
280     ForceFieldOptions& fopts = forceField_->getForceFieldOptions();
281 gezelter 1126
282 gezelter 1576 // Force fields can set options on how to scale van der Waals and electrostatic
283     // interactions for atoms connected via bonds, bends and torsions
284     // in this case the topological distance between atoms is:
285     // 0 = topologically unconnected
286     // 1 = bonded together
287     // 2 = connected via a bend
288     // 3 = connected via a torsion
289    
290     vdwScale_.reserve(4);
291     fill(vdwScale_.begin(), vdwScale_.end(), 0.0);
292    
293     electrostaticScale_.reserve(4);
294     fill(electrostaticScale_.begin(), electrostaticScale_.end(), 0.0);
295    
296     vdwScale_[0] = 1.0;
297     vdwScale_[1] = fopts.getvdw12scale();
298     vdwScale_[2] = fopts.getvdw13scale();
299     vdwScale_[3] = fopts.getvdw14scale();
300    
301     electrostaticScale_[0] = 1.0;
302     electrostaticScale_[1] = fopts.getelectrostatic12scale();
303     electrostaticScale_[2] = fopts.getelectrostatic13scale();
304     electrostaticScale_[3] = fopts.getelectrostatic14scale();
305    
306     fDecomp_->distributeInitialData();
307    
308     initialized_ = true;
309    
310     }
311    
312     void ForceManager::calcForces() {
313    
314     if (!initialized_) initialize();
315    
316 gezelter 1544 preCalculation();
317 gezelter 1546 shortRangeInteractions();
318     longRangeInteractions();
319 gezelter 1576 postCalculation();
320 gezelter 507 }
321 gezelter 1126
322 gezelter 507 void ForceManager::preCalculation() {
323 gezelter 246 SimInfo::MoleculeIterator mi;
324     Molecule* mol;
325     Molecule::AtomIterator ai;
326     Atom* atom;
327     Molecule::RigidBodyIterator rbIter;
328     RigidBody* rb;
329 gezelter 1540 Molecule::CutoffGroupIterator ci;
330     CutoffGroup* cg;
331 gezelter 246
332     // forces are zeroed here, before any are accumulated.
333 chuckv 1245
334 gezelter 1126 for (mol = info_->beginMolecule(mi); mol != NULL;
335     mol = info_->nextMolecule(mi)) {
336 gezelter 507 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
337     atom->zeroForcesAndTorques();
338     }
339 chuckv 1245
340 gezelter 507 //change the positions of atoms which belong to the rigidbodies
341 gezelter 1126 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
342     rb = mol->nextRigidBody(rbIter)) {
343 gezelter 507 rb->zeroForcesAndTorques();
344     }
345 gezelter 1540
346     if(info_->getNGlobalCutoffGroups() != info_->getNGlobalAtoms()){
347     for(cg = mol->beginCutoffGroup(ci); cg != NULL;
348     cg = mol->nextCutoffGroup(ci)) {
349     //calculate the center of mass of cutoff group
350     cg->updateCOM();
351     }
352     }
353 gezelter 246 }
354 gezelter 1540
355 gezelter 1126 // Zero out the stress tensor
356     tau *= 0.0;
357    
358 gezelter 507 }
359 gezelter 1126
360 gezelter 1546 void ForceManager::shortRangeInteractions() {
361 gezelter 246 Molecule* mol;
362     RigidBody* rb;
363     Bond* bond;
364     Bend* bend;
365     Torsion* torsion;
366 cli2 1275 Inversion* inversion;
367 gezelter 246 SimInfo::MoleculeIterator mi;
368     Molecule::RigidBodyIterator rbIter;
369     Molecule::BondIterator bondIter;;
370     Molecule::BendIterator bendIter;
371     Molecule::TorsionIterator torsionIter;
372 cli2 1275 Molecule::InversionIterator inversionIter;
373 tim 963 RealType bondPotential = 0.0;
374     RealType bendPotential = 0.0;
375     RealType torsionPotential = 0.0;
376 cli2 1275 RealType inversionPotential = 0.0;
377 gezelter 246
378     //calculate short range interactions
379 gezelter 1126 for (mol = info_->beginMolecule(mi); mol != NULL;
380     mol = info_->nextMolecule(mi)) {
381 gezelter 246
382 gezelter 507 //change the positions of atoms which belong to the rigidbodies
383 gezelter 1126 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
384     rb = mol->nextRigidBody(rbIter)) {
385     rb->updateAtoms();
386 gezelter 507 }
387 gezelter 246
388 gezelter 1126 for (bond = mol->beginBond(bondIter); bond != NULL;
389     bond = mol->nextBond(bondIter)) {
390 tim 749 bond->calcForce();
391     bondPotential += bond->getPotential();
392 gezelter 507 }
393 gezelter 246
394 gezelter 1126 for (bend = mol->beginBend(bendIter); bend != NULL;
395     bend = mol->nextBend(bendIter)) {
396    
397     RealType angle;
398     bend->calcForce(angle);
399     RealType currBendPot = bend->getPotential();
400 gezelter 1448
401 gezelter 1126 bendPotential += bend->getPotential();
402 gezelter 1545 map<Bend*, BendDataSet>::iterator i = bendDataSets.find(bend);
403 gezelter 1126 if (i == bendDataSets.end()) {
404     BendDataSet dataSet;
405     dataSet.prev.angle = dataSet.curr.angle = angle;
406     dataSet.prev.potential = dataSet.curr.potential = currBendPot;
407     dataSet.deltaV = 0.0;
408 gezelter 1545 bendDataSets.insert(map<Bend*, BendDataSet>::value_type(bend, dataSet));
409 gezelter 1126 }else {
410     i->second.prev.angle = i->second.curr.angle;
411     i->second.prev.potential = i->second.curr.potential;
412     i->second.curr.angle = angle;
413     i->second.curr.potential = currBendPot;
414     i->second.deltaV = fabs(i->second.curr.potential -
415     i->second.prev.potential);
416     }
417 gezelter 507 }
418 gezelter 1126
419     for (torsion = mol->beginTorsion(torsionIter); torsion != NULL;
420     torsion = mol->nextTorsion(torsionIter)) {
421 tim 963 RealType angle;
422 gezelter 1126 torsion->calcForce(angle);
423 tim 963 RealType currTorsionPot = torsion->getPotential();
424 gezelter 1126 torsionPotential += torsion->getPotential();
425 gezelter 1545 map<Torsion*, TorsionDataSet>::iterator i = torsionDataSets.find(torsion);
426 gezelter 1126 if (i == torsionDataSets.end()) {
427     TorsionDataSet dataSet;
428     dataSet.prev.angle = dataSet.curr.angle = angle;
429     dataSet.prev.potential = dataSet.curr.potential = currTorsionPot;
430     dataSet.deltaV = 0.0;
431 gezelter 1545 torsionDataSets.insert(map<Torsion*, TorsionDataSet>::value_type(torsion, dataSet));
432 gezelter 1126 }else {
433     i->second.prev.angle = i->second.curr.angle;
434     i->second.prev.potential = i->second.curr.potential;
435     i->second.curr.angle = angle;
436     i->second.curr.potential = currTorsionPot;
437     i->second.deltaV = fabs(i->second.curr.potential -
438     i->second.prev.potential);
439     }
440     }
441 gezelter 1545
442 cli2 1275 for (inversion = mol->beginInversion(inversionIter);
443     inversion != NULL;
444     inversion = mol->nextInversion(inversionIter)) {
445     RealType angle;
446     inversion->calcForce(angle);
447     RealType currInversionPot = inversion->getPotential();
448     inversionPotential += inversion->getPotential();
449 gezelter 1545 map<Inversion*, InversionDataSet>::iterator i = inversionDataSets.find(inversion);
450 cli2 1275 if (i == inversionDataSets.end()) {
451     InversionDataSet dataSet;
452     dataSet.prev.angle = dataSet.curr.angle = angle;
453     dataSet.prev.potential = dataSet.curr.potential = currInversionPot;
454     dataSet.deltaV = 0.0;
455 gezelter 1545 inversionDataSets.insert(map<Inversion*, InversionDataSet>::value_type(inversion, dataSet));
456 cli2 1275 }else {
457     i->second.prev.angle = i->second.curr.angle;
458     i->second.prev.potential = i->second.curr.potential;
459     i->second.curr.angle = angle;
460     i->second.curr.potential = currInversionPot;
461     i->second.deltaV = fabs(i->second.curr.potential -
462     i->second.prev.potential);
463     }
464     }
465 gezelter 246 }
466    
467 gezelter 1126 RealType shortRangePotential = bondPotential + bendPotential +
468 cli2 1275 torsionPotential + inversionPotential;
469 gezelter 246 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
470     curSnapshot->statData[Stats::SHORT_RANGE_POTENTIAL] = shortRangePotential;
471 tim 665 curSnapshot->statData[Stats::BOND_POTENTIAL] = bondPotential;
472     curSnapshot->statData[Stats::BEND_POTENTIAL] = bendPotential;
473     curSnapshot->statData[Stats::DIHEDRAL_POTENTIAL] = torsionPotential;
474 gezelter 1545 curSnapshot->statData[Stats::INVERSION_POTENTIAL] = inversionPotential;
475 gezelter 507 }
476 gezelter 1126
477 gezelter 1546 void ForceManager::longRangeInteractions() {
478 gezelter 1581
479 gezelter 1545 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
480     DataStorage* config = &(curSnapshot->atomData);
481     DataStorage* cgConfig = &(curSnapshot->cgData);
482    
483 gezelter 1581 //calculate the center of mass of cutoff group
484    
485     SimInfo::MoleculeIterator mi;
486     Molecule* mol;
487     Molecule::CutoffGroupIterator ci;
488     CutoffGroup* cg;
489    
490     if(info_->getNCutoffGroups() > 0){
491     for (mol = info_->beginMolecule(mi); mol != NULL;
492     mol = info_->nextMolecule(mi)) {
493     for(cg = mol->beginCutoffGroup(ci); cg != NULL;
494     cg = mol->nextCutoffGroup(ci)) {
495     cg->updateCOM();
496     }
497     }
498     } else {
499     // center of mass of the group is the same as position of the atom
500     // if cutoff group does not exist
501     cgConfig->position = config->position;
502     }
503    
504 gezelter 1575 fDecomp_->zeroWorkArrays();
505 gezelter 1549 fDecomp_->distributeData();
506 gezelter 1579
507     int cg1, cg2, atom1, atom2, topoDist;
508     Vector3d d_grp, dag, d;
509     RealType rgrpsq, rgrp, r2, r;
510     RealType electroMult, vdwMult;
511 gezelter 1549 RealType vij;
512 gezelter 1581 Vector3d fij, fg, f1;
513 gezelter 1576 tuple3<RealType, RealType, RealType> cuts;
514 gezelter 1545 RealType rCutSq;
515     bool in_switching_region;
516     RealType sw, dswdr, swderiv;
517 gezelter 1549 vector<int> atomListColumn, atomListRow, atomListLocal;
518 gezelter 1545 InteractionData idat;
519 gezelter 1546 SelfData sdat;
520     RealType mf;
521 gezelter 1575 potVec pot(0.0);
522     potVec longRangePotential(0.0);
523     RealType lrPot;
524 gezelter 1579 RealType vpair;
525 gezelter 1544
526 gezelter 1545 int loopStart, loopEnd;
527 gezelter 1544
528 gezelter 1581 idat.vdwMult = &vdwMult;
529     idat.electroMult = &electroMult;
530     idat.pot = &pot;
531     idat.vpair = &vpair;
532     idat.f1 = &f1;
533     idat.sw = &sw;
534    
535 gezelter 1545 loopEnd = PAIR_LOOP;
536 gezelter 1546 if (info_->requiresPrepair() ) {
537 gezelter 1545 loopStart = PREPAIR_LOOP;
538     } else {
539     loopStart = PAIR_LOOP;
540     }
541 gezelter 1579
542     for (int iLoop = loopStart; iLoop <= loopEnd; iLoop++) {
543    
544 gezelter 1545 if (iLoop == loopStart) {
545 gezelter 1549 bool update_nlist = fDecomp_->checkNeighborList();
546 gezelter 1545 if (update_nlist)
547 gezelter 1549 neighborList = fDecomp_->buildNeighborList();
548 gezelter 1579 }
549    
550 gezelter 1545 for (vector<pair<int, int> >::iterator it = neighborList.begin();
551     it != neighborList.end(); ++it) {
552 gezelter 1579
553 gezelter 1545 cg1 = (*it).first;
554     cg2 = (*it).second;
555 gezelter 1576
556     cuts = fDecomp_->getGroupCutoffs(cg1, cg2);
557 gezelter 1545
558 gezelter 1549 d_grp = fDecomp_->getIntergroupVector(cg1, cg2);
559 gezelter 1545 curSnapshot->wrapVector(d_grp);
560     rgrpsq = d_grp.lengthSquare();
561    
562 gezelter 1576 rCutSq = cuts.second;
563    
564 gezelter 1545 if (rgrpsq < rCutSq) {
565 gezelter 1579 idat.rcut = &cuts.first;
566 gezelter 1545 if (iLoop == PAIR_LOOP) {
567 gezelter 1546 vij *= 0.0;
568 gezelter 1545 fij = V3Zero;
569     }
570    
571 gezelter 1579 in_switching_region = switcher_->getSwitch(rgrpsq, sw, dswdr,
572 gezelter 1576 rgrp);
573    
574 gezelter 1549 atomListRow = fDecomp_->getAtomsInGroupRow(cg1);
575     atomListColumn = fDecomp_->getAtomsInGroupColumn(cg2);
576 gezelter 1545
577 gezelter 1549 for (vector<int>::iterator ia = atomListRow.begin();
578     ia != atomListRow.end(); ++ia) {
579 gezelter 1545 atom1 = (*ia);
580    
581 gezelter 1549 for (vector<int>::iterator jb = atomListColumn.begin();
582     jb != atomListColumn.end(); ++jb) {
583 gezelter 1545 atom2 = (*jb);
584    
585 gezelter 1549 if (!fDecomp_->skipAtomPair(atom1, atom2)) {
586 gezelter 1545
587 gezelter 1579 vpair = 0.0;
588 gezelter 1581 f1 = V3Zero;
589 gezelter 1575
590 gezelter 1581 fDecomp_->fillInteractionData(idat, atom1, atom2);
591 gezelter 1579
592     topoDist = fDecomp_->getTopologicalDistance(atom1, atom2);
593     vdwMult = vdwScale_[topoDist];
594     electroMult = electrostaticScale_[topoDist];
595 gezelter 1546
596 gezelter 1549 if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
597 gezelter 1579 idat.d = &d_grp;
598     idat.r2 = &rgrpsq;
599 gezelter 1545 } else {
600 gezelter 1579 d = fDecomp_->getInteratomicVector(atom1, atom2);
601     curSnapshot->wrapVector( d );
602     r2 = d.lengthSquare();
603     idat.d = &d;
604     idat.r2 = &r2;
605 gezelter 1545 }
606    
607 gezelter 1581 r = sqrt( *(idat.r2) );
608 gezelter 1579 idat.rij = &r;
609 gezelter 1546
610 gezelter 1545 if (iLoop == PREPAIR_LOOP) {
611     interactionMan_->doPrePair(idat);
612     } else {
613     interactionMan_->doPair(idat);
614 gezelter 1575 fDecomp_->unpackInteractionData(idat, atom1, atom2);
615 gezelter 1581 vij += vpair;
616     fij += f1;
617     tau -= outProduct( *(idat.d), f1);
618 gezelter 1545 }
619     }
620     }
621     }
622    
623     if (iLoop == PAIR_LOOP) {
624     if (in_switching_region) {
625     swderiv = vij * dswdr / rgrp;
626     fg = swderiv * d_grp;
627    
628     fij += fg;
629    
630 gezelter 1549 if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
631 gezelter 1554 tau -= outProduct( *(idat.d), fg);
632 gezelter 1545 }
633    
634 gezelter 1549 for (vector<int>::iterator ia = atomListRow.begin();
635     ia != atomListRow.end(); ++ia) {
636 gezelter 1545 atom1 = (*ia);
637 gezelter 1569 mf = fDecomp_->getMassFactorRow(atom1);
638 gezelter 1545 // fg is the force on atom ia due to cutoff group's
639     // presence in switching region
640     fg = swderiv * d_grp * mf;
641 gezelter 1549 fDecomp_->addForceToAtomRow(atom1, fg);
642 gezelter 1545
643 gezelter 1549 if (atomListRow.size() > 1) {
644 gezelter 1546 if (info_->usesAtomicVirial()) {
645 gezelter 1545 // find the distance between the atom
646     // and the center of the cutoff group:
647 gezelter 1549 dag = fDecomp_->getAtomToGroupVectorRow(atom1, cg1);
648 gezelter 1545 tau -= outProduct(dag, fg);
649     }
650     }
651     }
652 gezelter 1549 for (vector<int>::iterator jb = atomListColumn.begin();
653     jb != atomListColumn.end(); ++jb) {
654 gezelter 1545 atom2 = (*jb);
655 gezelter 1569 mf = fDecomp_->getMassFactorColumn(atom2);
656 gezelter 1545 // fg is the force on atom jb due to cutoff group's
657     // presence in switching region
658     fg = -swderiv * d_grp * mf;
659 gezelter 1549 fDecomp_->addForceToAtomColumn(atom2, fg);
660 gezelter 1545
661 gezelter 1549 if (atomListColumn.size() > 1) {
662 gezelter 1546 if (info_->usesAtomicVirial()) {
663 gezelter 1545 // find the distance between the atom
664     // and the center of the cutoff group:
665 gezelter 1549 dag = fDecomp_->getAtomToGroupVectorColumn(atom2, cg2);
666 gezelter 1545 tau -= outProduct(dag, fg);
667     }
668     }
669     }
670     }
671     //if (!SIM_uses_AtomicVirial) {
672     // tau -= outProduct(d_grp, fij);
673     //}
674     }
675     }
676     }
677    
678     if (iLoop == PREPAIR_LOOP) {
679 gezelter 1546 if (info_->requiresPrepair()) {
680 gezelter 1549 fDecomp_->collectIntermediateData();
681 gezelter 1570
682     for (int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
683 gezelter 1581 fDecomp_->fillSelfData(sdat, atom1);
684 gezelter 1545 interactionMan_->doPreForce(sdat);
685     }
686 gezelter 1570
687 gezelter 1549 fDecomp_->distributeIntermediateData();
688 gezelter 1545 }
689     }
690    
691 gezelter 1544 }
692 gezelter 1545
693 gezelter 1549 fDecomp_->collectData();
694 gezelter 1545
695 gezelter 1570 if ( info_->requiresSkipCorrection() ) {
696    
697     for (int atom1 = 0; atom1 < fDecomp_->getNAtomsInRow(); atom1++) {
698 gezelter 1544
699 gezelter 1579 vector<int> skipList = fDecomp_->getSkipsForAtom( atom1 );
700 gezelter 1570
701     for (vector<int>::iterator jb = skipList.begin();
702     jb != skipList.end(); ++jb) {
703    
704     atom2 = (*jb);
705 gezelter 1581 fDecomp_->fillSkipData(idat, atom1, atom2);
706 gezelter 1570 interactionMan_->doSkipCorrection(idat);
707    
708 gezelter 1545 }
709     }
710 gezelter 246 }
711 gezelter 1570
712     if (info_->requiresSelfCorrection()) {
713 gezelter 1545
714 gezelter 1570 for (int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
715 gezelter 1581 fDecomp_->fillSelfData(sdat, atom1);
716 gezelter 1570 interactionMan_->doSelfCorrection(sdat);
717     }
718    
719     }
720    
721 gezelter 1575 longRangePotential = fDecomp_->getLongRangePotential();
722     lrPot = longRangePotential.sum();
723    
724 gezelter 246 //store the tau and long range potential
725 chuckv 664 curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL] = lrPot;
726 gezelter 1550 curSnapshot->statData[Stats::VANDERWAALS_POTENTIAL] = longRangePotential[VANDERWAALS_FAMILY];
727     curSnapshot->statData[Stats::ELECTROSTATIC_POTENTIAL] = longRangePotential[ELECTROSTATIC_FAMILY];
728 gezelter 507 }
729 gezelter 246
730 gezelter 1126
731 gezelter 1464 void ForceManager::postCalculation() {
732 gezelter 246 SimInfo::MoleculeIterator mi;
733     Molecule* mol;
734     Molecule::RigidBodyIterator rbIter;
735     RigidBody* rb;
736 gezelter 1126 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
737 gezelter 246
738     // collect the atomic forces onto rigid bodies
739 gezelter 1126
740     for (mol = info_->beginMolecule(mi); mol != NULL;
741     mol = info_->nextMolecule(mi)) {
742     for (rb = mol->beginRigidBody(rbIter); rb != NULL;
743     rb = mol->nextRigidBody(rbIter)) {
744 gezelter 1464 Mat3x3d rbTau = rb->calcForcesAndTorquesAndVirial();
745     tau += rbTau;
746 gezelter 507 }
747 gezelter 1126 }
748 gezelter 1464
749 gezelter 1126 #ifdef IS_MPI
750 gezelter 1464 Mat3x3d tmpTau(tau);
751     MPI_Allreduce(tmpTau.getArrayPointer(), tau.getArrayPointer(),
752     9, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
753 gezelter 1126 #endif
754 gezelter 1464 curSnapshot->statData.setTau(tau);
755 gezelter 507 }
756 gezelter 246
757 gezelter 1390 } //end namespace OpenMD

Properties

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