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

Properties

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