ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/ForceManager.cpp
Revision: 1711
Committed: Sat May 19 02:58:35 2012 UTC (12 years, 11 months ago) by gezelter
File size: 31759 byte(s)
Log Message:
Some fixes for DataStorage issues.  Removed outdated zangle stuff that
has been replaced by the more modern restraints.

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

Properties

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