ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/ForceManager.cpp
Revision: 1847
Committed: Thu Jan 31 17:57:07 2013 UTC (12 years, 3 months ago) by gezelter
File size: 35165 byte(s)
Log Message:
Fixed a field non-zeroing bug

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     * @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 jmarr 1780 #include "perturbations/ElectricField.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 jmichalk 1754 rCut_ = 0.0; //Needs a value for a later max() call;
114 gezelter 1755
115 gezelter 1613 if (simParams_->haveMDfileVersion())
116     mdFileVersion = simParams_->getMDfileVersion();
117     else
118     mdFileVersion = 0;
119    
120 gezelter 1576 if (simParams_->haveCutoffRadius()) {
121     rCut_ = simParams_->getCutoffRadius();
122     } else {
123     if (info_->usesElectrostaticAtoms()) {
124     sprintf(painCave.errMsg,
125     "ForceManager::setupCutoffs: No value was set for the cutoffRadius.\n"
126     "\tOpenMD will use a default value of 12.0 angstroms"
127     "\tfor the cutoffRadius.\n");
128     painCave.isFatal = 0;
129     painCave.severity = OPENMD_INFO;
130     simError();
131     rCut_ = 12.0;
132     } else {
133     RealType thisCut;
134     set<AtomType*>::iterator i;
135     set<AtomType*> atomTypes;
136     atomTypes = info_->getSimulatedAtomTypes();
137     for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
138     thisCut = interactionMan_->getSuggestedCutoffRadius((*i));
139     rCut_ = max(thisCut, rCut_);
140     }
141     sprintf(painCave.errMsg,
142     "ForceManager::setupCutoffs: No value was set for the cutoffRadius.\n"
143     "\tOpenMD will use %lf angstroms.\n",
144     rCut_);
145     painCave.isFatal = 0;
146     painCave.severity = OPENMD_INFO;
147     simError();
148 gezelter 1579 }
149 gezelter 1576 }
150    
151 gezelter 1583 fDecomp_->setUserCutoff(rCut_);
152 gezelter 1584 interactionMan_->setCutoffRadius(rCut_);
153 gezelter 1583
154 gezelter 1576 map<string, CutoffMethod> stringToCutoffMethod;
155     stringToCutoffMethod["HARD"] = HARD;
156     stringToCutoffMethod["SWITCHED"] = SWITCHED;
157     stringToCutoffMethod["SHIFTED_POTENTIAL"] = SHIFTED_POTENTIAL;
158     stringToCutoffMethod["SHIFTED_FORCE"] = SHIFTED_FORCE;
159 gezelter 1545
160 gezelter 1576 if (simParams_->haveCutoffMethod()) {
161     string cutMeth = toUpperCopy(simParams_->getCutoffMethod());
162     map<string, CutoffMethod>::iterator i;
163     i = stringToCutoffMethod.find(cutMeth);
164     if (i == stringToCutoffMethod.end()) {
165     sprintf(painCave.errMsg,
166     "ForceManager::setupCutoffs: Could not find chosen cutoffMethod %s\n"
167     "\tShould be one of: "
168     "HARD, SWITCHED, SHIFTED_POTENTIAL, or SHIFTED_FORCE\n",
169     cutMeth.c_str());
170     painCave.isFatal = 1;
171     painCave.severity = OPENMD_ERROR;
172     simError();
173     } else {
174     cutoffMethod_ = i->second;
175     }
176     } else {
177 gezelter 1616 if (mdFileVersion > 1) {
178     sprintf(painCave.errMsg,
179     "ForceManager::setupCutoffs: No value was set for the cutoffMethod.\n"
180     "\tOpenMD will use SHIFTED_FORCE.\n");
181     painCave.isFatal = 0;
182     painCave.severity = OPENMD_INFO;
183     simError();
184     cutoffMethod_ = SHIFTED_FORCE;
185     } else {
186     // handle the case where the old file version was in play
187     // (there should be no cutoffMethod, so we have to deduce it
188     // from other data).
189    
190     sprintf(painCave.errMsg,
191     "ForceManager::setupCutoffs : DEPRECATED FILE FORMAT!\n"
192     "\tOpenMD found a file which does not set a cutoffMethod.\n"
193     "\tOpenMD will attempt to deduce a cutoffMethod using the\n"
194     "\tbehavior of the older (version 1) code. To remove this\n"
195     "\twarning, add an explicit cutoffMethod and change the top\n"
196     "\tof the file so that it begins with <OpenMD version=2>\n");
197     painCave.isFatal = 0;
198     painCave.severity = OPENMD_WARNING;
199     simError();
200    
201     // The old file version tethered the shifting behavior to the
202     // electrostaticSummationMethod keyword.
203    
204     if (simParams_->haveElectrostaticSummationMethod()) {
205 gezelter 1710 string myMethod = simParams_->getElectrostaticSummationMethod();
206 gezelter 1616 toUpper(myMethod);
207    
208     if (myMethod == "SHIFTED_POTENTIAL") {
209     cutoffMethod_ = SHIFTED_POTENTIAL;
210     } else if (myMethod == "SHIFTED_FORCE") {
211     cutoffMethod_ = SHIFTED_FORCE;
212     }
213    
214     if (simParams_->haveSwitchingRadius())
215     rSwitch_ = simParams_->getSwitchingRadius();
216    
217     if (myMethod == "SHIFTED_POTENTIAL" || myMethod == "SHIFTED_FORCE") {
218     if (simParams_->haveSwitchingRadius()){
219     sprintf(painCave.errMsg,
220     "ForceManager::setupCutoffs : DEPRECATED ERROR MESSAGE\n"
221     "\tA value was set for the switchingRadius\n"
222     "\teven though the electrostaticSummationMethod was\n"
223     "\tset to %s\n", myMethod.c_str());
224     painCave.severity = OPENMD_WARNING;
225     painCave.isFatal = 1;
226     simError();
227     }
228     }
229     if (abs(rCut_ - rSwitch_) < 0.0001) {
230     if (cutoffMethod_ == SHIFTED_FORCE) {
231     sprintf(painCave.errMsg,
232     "ForceManager::setupCutoffs : DEPRECATED BEHAVIOR\n"
233     "\tcutoffRadius and switchingRadius are set to the\n"
234     "\tsame value. OpenMD will use shifted force\n"
235     "\tpotentials instead of switching functions.\n");
236     painCave.isFatal = 0;
237     painCave.severity = OPENMD_WARNING;
238     simError();
239     } else {
240     cutoffMethod_ = SHIFTED_POTENTIAL;
241     sprintf(painCave.errMsg,
242     "ForceManager::setupCutoffs : DEPRECATED BEHAVIOR\n"
243     "\tcutoffRadius and switchingRadius are set to the\n"
244     "\tsame value. OpenMD will use shifted potentials\n"
245     "\tinstead of switching functions.\n");
246     painCave.isFatal = 0;
247     painCave.severity = OPENMD_WARNING;
248     simError();
249     }
250     }
251     }
252     }
253 gezelter 1576 }
254    
255     map<string, CutoffPolicy> stringToCutoffPolicy;
256     stringToCutoffPolicy["MIX"] = MIX;
257     stringToCutoffPolicy["MAX"] = MAX;
258     stringToCutoffPolicy["TRADITIONAL"] = TRADITIONAL;
259    
260 gezelter 1710 string cutPolicy;
261 gezelter 1576 if (forceFieldOptions_.haveCutoffPolicy()){
262     cutPolicy = forceFieldOptions_.getCutoffPolicy();
263     }else if (simParams_->haveCutoffPolicy()) {
264     cutPolicy = simParams_->getCutoffPolicy();
265     }
266    
267     if (!cutPolicy.empty()){
268     toUpper(cutPolicy);
269     map<string, CutoffPolicy>::iterator i;
270     i = stringToCutoffPolicy.find(cutPolicy);
271    
272     if (i == stringToCutoffPolicy.end()) {
273     sprintf(painCave.errMsg,
274     "ForceManager::setupCutoffs: Could not find chosen cutoffPolicy %s\n"
275     "\tShould be one of: "
276     "MIX, MAX, or TRADITIONAL\n",
277     cutPolicy.c_str());
278     painCave.isFatal = 1;
279     painCave.severity = OPENMD_ERROR;
280     simError();
281     } else {
282     cutoffPolicy_ = i->second;
283     }
284     } else {
285     sprintf(painCave.errMsg,
286     "ForceManager::setupCutoffs: No value was set for the cutoffPolicy.\n"
287     "\tOpenMD will use TRADITIONAL.\n");
288     painCave.isFatal = 0;
289     painCave.severity = OPENMD_INFO;
290     simError();
291     cutoffPolicy_ = TRADITIONAL;
292     }
293 gezelter 1587
294 gezelter 1579 fDecomp_->setCutoffPolicy(cutoffPolicy_);
295 gezelter 1587
296     // create the switching function object:
297 gezelter 1576
298 gezelter 1577 switcher_ = new SwitchingFunction();
299 gezelter 1587
300     if (cutoffMethod_ == SWITCHED) {
301     if (simParams_->haveSwitchingRadius()) {
302     rSwitch_ = simParams_->getSwitchingRadius();
303     if (rSwitch_ > rCut_) {
304     sprintf(painCave.errMsg,
305     "ForceManager::setupCutoffs: switchingRadius (%f) is larger "
306     "than the cutoffRadius(%f)\n", rSwitch_, rCut_);
307     painCave.isFatal = 1;
308     painCave.severity = OPENMD_ERROR;
309     simError();
310     }
311     } else {
312     rSwitch_ = 0.85 * rCut_;
313 gezelter 1576 sprintf(painCave.errMsg,
314 gezelter 1587 "ForceManager::setupCutoffs: No value was set for the switchingRadius.\n"
315     "\tOpenMD will use a default value of 85 percent of the cutoffRadius.\n"
316     "\tswitchingRadius = %f. for this simulation\n", rSwitch_);
317     painCave.isFatal = 0;
318     painCave.severity = OPENMD_WARNING;
319 gezelter 1576 simError();
320     }
321 gezelter 1587 } else {
322 gezelter 1618 if (mdFileVersion > 1) {
323     // throw an error if we define a switching radius and don't need one.
324     // older file versions should not do this.
325     if (simParams_->haveSwitchingRadius()) {
326     map<string, CutoffMethod>::const_iterator it;
327     string theMeth;
328     for (it = stringToCutoffMethod.begin();
329     it != stringToCutoffMethod.end(); ++it) {
330     if (it->second == cutoffMethod_) {
331     theMeth = it->first;
332     break;
333     }
334 gezelter 1587 }
335 gezelter 1618 sprintf(painCave.errMsg,
336     "ForceManager::setupCutoffs: the cutoffMethod (%s)\n"
337     "\tis not set to SWITCHED, so switchingRadius value\n"
338     "\twill be ignored for this simulation\n", theMeth.c_str());
339     painCave.isFatal = 0;
340     painCave.severity = OPENMD_WARNING;
341     simError();
342 gezelter 1587 }
343     }
344     rSwitch_ = rCut_;
345     }
346 gezelter 1576
347 gezelter 1577 // Default to cubic switching function.
348     sft_ = cubic;
349 gezelter 1576 if (simParams_->haveSwitchingFunctionType()) {
350     string funcType = simParams_->getSwitchingFunctionType();
351     toUpper(funcType);
352     if (funcType == "CUBIC") {
353     sft_ = cubic;
354     } else {
355     if (funcType == "FIFTH_ORDER_POLYNOMIAL") {
356     sft_ = fifth_order_poly;
357     } else {
358     // throw error
359     sprintf( painCave.errMsg,
360     "ForceManager::setupSwitching : Unknown switchingFunctionType. (Input file specified %s .)\n"
361     "\tswitchingFunctionType must be one of: "
362     "\"cubic\" or \"fifth_order_polynomial\".",
363     funcType.c_str() );
364     painCave.isFatal = 1;
365     painCave.severity = OPENMD_ERROR;
366     simError();
367     }
368     }
369     }
370     switcher_->setSwitchType(sft_);
371     switcher_->setSwitch(rSwitch_, rCut_);
372     }
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 gezelter 1723 doHeatFlux_ = info_->getSimParams()->getPrintHeatFlux();
394     if (doHeatFlux_) doParticlePot_ = true;
395 gezelter 1787
396     doElectricField_ = info_->getSimParams()->getOutputElectricField();
397 gezelter 1711
398 gezelter 246 }
399 gezelter 1576
400     ForceFieldOptions& fopts = forceField_->getForceFieldOptions();
401 gezelter 1126
402 gezelter 1590 // Force fields can set options on how to scale van der Waals and
403     // electrostatic interactions for atoms connected via bonds, bends
404     // and torsions in this case the topological distance between
405     // atoms is:
406 gezelter 1576 // 0 = topologically unconnected
407     // 1 = bonded together
408     // 2 = connected via a bend
409     // 3 = connected via a torsion
410    
411     vdwScale_.reserve(4);
412     fill(vdwScale_.begin(), vdwScale_.end(), 0.0);
413    
414     electrostaticScale_.reserve(4);
415     fill(electrostaticScale_.begin(), electrostaticScale_.end(), 0.0);
416    
417     vdwScale_[0] = 1.0;
418     vdwScale_[1] = fopts.getvdw12scale();
419     vdwScale_[2] = fopts.getvdw13scale();
420     vdwScale_[3] = fopts.getvdw14scale();
421    
422     electrostaticScale_[0] = 1.0;
423     electrostaticScale_[1] = fopts.getelectrostatic12scale();
424     electrostaticScale_[2] = fopts.getelectrostatic13scale();
425     electrostaticScale_[3] = fopts.getelectrostatic14scale();
426    
427 jmarr 1780 if (info_->getSimParams()->haveElectricField()) {
428     ElectricField* eField = new ElectricField(info_);
429     perturbations_.push_back(eField);
430     }
431    
432 gezelter 1576 fDecomp_->distributeInitialData();
433    
434     initialized_ = true;
435    
436     }
437    
438     void ForceManager::calcForces() {
439    
440     if (!initialized_) initialize();
441    
442 gezelter 1544 preCalculation();
443 gezelter 1546 shortRangeInteractions();
444     longRangeInteractions();
445 gezelter 1576 postCalculation();
446 gezelter 507 }
447 gezelter 1126
448 gezelter 507 void ForceManager::preCalculation() {
449 gezelter 246 SimInfo::MoleculeIterator mi;
450     Molecule* mol;
451     Molecule::AtomIterator ai;
452     Atom* atom;
453     Molecule::RigidBodyIterator rbIter;
454     RigidBody* rb;
455 gezelter 1540 Molecule::CutoffGroupIterator ci;
456     CutoffGroup* cg;
457 gezelter 246
458 gezelter 1764 // forces and potentials are zeroed here, before any are
459     // accumulated.
460 chuckv 1245
461 gezelter 1764 Snapshot* snap = info_->getSnapshotManager()->getCurrentSnapshot();
462    
463     snap->setBondPotential(0.0);
464     snap->setBendPotential(0.0);
465     snap->setTorsionPotential(0.0);
466     snap->setInversionPotential(0.0);
467    
468     potVec zeroPot(0.0);
469     snap->setLongRangePotential(zeroPot);
470     snap->setExcludedPotentials(zeroPot);
471    
472     snap->setRestraintPotential(0.0);
473     snap->setRawPotential(0.0);
474    
475 gezelter 1126 for (mol = info_->beginMolecule(mi); mol != NULL;
476     mol = info_->nextMolecule(mi)) {
477 gezelter 1590 for(atom = mol->beginAtom(ai); atom != NULL;
478     atom = mol->nextAtom(ai)) {
479 gezelter 507 atom->zeroForcesAndTorques();
480     }
481 gezelter 1590
482 gezelter 507 //change the positions of atoms which belong to the rigidbodies
483 gezelter 1126 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
484     rb = mol->nextRigidBody(rbIter)) {
485 gezelter 507 rb->zeroForcesAndTorques();
486     }
487 gezelter 1590
488 gezelter 1540 if(info_->getNGlobalCutoffGroups() != info_->getNGlobalAtoms()){
489     for(cg = mol->beginCutoffGroup(ci); cg != NULL;
490     cg = mol->nextCutoffGroup(ci)) {
491     //calculate the center of mass of cutoff group
492     cg->updateCOM();
493     }
494     }
495 gezelter 246 }
496 gezelter 1590
497 gezelter 1126 // Zero out the stress tensor
498 gezelter 1723 stressTensor *= 0.0;
499     // Zero out the heatFlux
500 gezelter 1744 fDecomp_->setHeatFlux( Vector3d(0.0) );
501 gezelter 507 }
502 gezelter 1126
503 gezelter 1546 void ForceManager::shortRangeInteractions() {
504 gezelter 246 Molecule* mol;
505     RigidBody* rb;
506     Bond* bond;
507     Bend* bend;
508     Torsion* torsion;
509 cli2 1275 Inversion* inversion;
510 gezelter 246 SimInfo::MoleculeIterator mi;
511     Molecule::RigidBodyIterator rbIter;
512     Molecule::BondIterator bondIter;;
513     Molecule::BendIterator bendIter;
514     Molecule::TorsionIterator torsionIter;
515 cli2 1275 Molecule::InversionIterator inversionIter;
516 tim 963 RealType bondPotential = 0.0;
517     RealType bendPotential = 0.0;
518     RealType torsionPotential = 0.0;
519 cli2 1275 RealType inversionPotential = 0.0;
520 gezelter 246
521     //calculate short range interactions
522 gezelter 1126 for (mol = info_->beginMolecule(mi); mol != NULL;
523     mol = info_->nextMolecule(mi)) {
524 gezelter 246
525 gezelter 507 //change the positions of atoms which belong to the rigidbodies
526 gezelter 1126 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
527     rb = mol->nextRigidBody(rbIter)) {
528     rb->updateAtoms();
529 gezelter 507 }
530 gezelter 246
531 gezelter 1126 for (bond = mol->beginBond(bondIter); bond != NULL;
532     bond = mol->nextBond(bondIter)) {
533 gezelter 1712 bond->calcForce(doParticlePot_);
534 tim 749 bondPotential += bond->getPotential();
535 gezelter 507 }
536 gezelter 246
537 gezelter 1126 for (bend = mol->beginBend(bendIter); bend != NULL;
538     bend = mol->nextBend(bendIter)) {
539    
540     RealType angle;
541 gezelter 1712 bend->calcForce(angle, doParticlePot_);
542 gezelter 1126 RealType currBendPot = bend->getPotential();
543 gezelter 1448
544 gezelter 1126 bendPotential += bend->getPotential();
545 gezelter 1545 map<Bend*, BendDataSet>::iterator i = bendDataSets.find(bend);
546 gezelter 1126 if (i == bendDataSets.end()) {
547     BendDataSet dataSet;
548     dataSet.prev.angle = dataSet.curr.angle = angle;
549     dataSet.prev.potential = dataSet.curr.potential = currBendPot;
550     dataSet.deltaV = 0.0;
551 gezelter 1590 bendDataSets.insert(map<Bend*, BendDataSet>::value_type(bend,
552     dataSet));
553 gezelter 1126 }else {
554     i->second.prev.angle = i->second.curr.angle;
555     i->second.prev.potential = i->second.curr.potential;
556     i->second.curr.angle = angle;
557     i->second.curr.potential = currBendPot;
558     i->second.deltaV = fabs(i->second.curr.potential -
559     i->second.prev.potential);
560     }
561 gezelter 507 }
562 gezelter 1126
563     for (torsion = mol->beginTorsion(torsionIter); torsion != NULL;
564     torsion = mol->nextTorsion(torsionIter)) {
565 tim 963 RealType angle;
566 gezelter 1712 torsion->calcForce(angle, doParticlePot_);
567 tim 963 RealType currTorsionPot = torsion->getPotential();
568 gezelter 1126 torsionPotential += torsion->getPotential();
569 gezelter 1545 map<Torsion*, TorsionDataSet>::iterator i = torsionDataSets.find(torsion);
570 gezelter 1126 if (i == torsionDataSets.end()) {
571     TorsionDataSet dataSet;
572     dataSet.prev.angle = dataSet.curr.angle = angle;
573     dataSet.prev.potential = dataSet.curr.potential = currTorsionPot;
574     dataSet.deltaV = 0.0;
575 gezelter 1545 torsionDataSets.insert(map<Torsion*, TorsionDataSet>::value_type(torsion, dataSet));
576 gezelter 1126 }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 = currTorsionPot;
581     i->second.deltaV = fabs(i->second.curr.potential -
582     i->second.prev.potential);
583     }
584     }
585 gezelter 1545
586 cli2 1275 for (inversion = mol->beginInversion(inversionIter);
587     inversion != NULL;
588     inversion = mol->nextInversion(inversionIter)) {
589     RealType angle;
590 gezelter 1712 inversion->calcForce(angle, doParticlePot_);
591 cli2 1275 RealType currInversionPot = inversion->getPotential();
592     inversionPotential += inversion->getPotential();
593 gezelter 1545 map<Inversion*, InversionDataSet>::iterator i = inversionDataSets.find(inversion);
594 cli2 1275 if (i == inversionDataSets.end()) {
595     InversionDataSet dataSet;
596     dataSet.prev.angle = dataSet.curr.angle = angle;
597     dataSet.prev.potential = dataSet.curr.potential = currInversionPot;
598     dataSet.deltaV = 0.0;
599 gezelter 1545 inversionDataSets.insert(map<Inversion*, InversionDataSet>::value_type(inversion, dataSet));
600 cli2 1275 }else {
601     i->second.prev.angle = i->second.curr.angle;
602     i->second.prev.potential = i->second.curr.potential;
603     i->second.curr.angle = angle;
604     i->second.curr.potential = currInversionPot;
605     i->second.deltaV = fabs(i->second.curr.potential -
606     i->second.prev.potential);
607     }
608     }
609 gezelter 246 }
610 gezelter 1760
611     #ifdef IS_MPI
612     // Collect from all nodes. This should eventually be moved into a
613     // SystemDecomposition, but this is a better place than in
614     // Thermo to do the collection.
615     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &bondPotential, 1, MPI::REALTYPE,
616     MPI::SUM);
617     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &bendPotential, 1, MPI::REALTYPE,
618     MPI::SUM);
619     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &torsionPotential, 1,
620     MPI::REALTYPE, MPI::SUM);
621     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &inversionPotential, 1,
622     MPI::REALTYPE, MPI::SUM);
623     #endif
624    
625     Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
626    
627     curSnapshot->setBondPotential(bondPotential);
628     curSnapshot->setBendPotential(bendPotential);
629     curSnapshot->setTorsionPotential(torsionPotential);
630     curSnapshot->setInversionPotential(inversionPotential);
631 gezelter 246
632 gezelter 1764 // RealType shortRangePotential = bondPotential + bendPotential +
633     // torsionPotential + inversionPotential;
634 gezelter 1760
635 gezelter 1764 // curSnapshot->setShortRangePotential(shortRangePotential);
636 gezelter 507 }
637 gezelter 1126
638 gezelter 1546 void ForceManager::longRangeInteractions() {
639 gezelter 1581
640 gezelter 1723
641 gezelter 1545 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
642     DataStorage* config = &(curSnapshot->atomData);
643     DataStorage* cgConfig = &(curSnapshot->cgData);
644    
645 gezelter 1581 //calculate the center of mass of cutoff group
646    
647     SimInfo::MoleculeIterator mi;
648     Molecule* mol;
649     Molecule::CutoffGroupIterator ci;
650     CutoffGroup* cg;
651    
652     if(info_->getNCutoffGroups() > 0){
653     for (mol = info_->beginMolecule(mi); mol != NULL;
654     mol = info_->nextMolecule(mi)) {
655     for(cg = mol->beginCutoffGroup(ci); cg != NULL;
656     cg = mol->nextCutoffGroup(ci)) {
657     cg->updateCOM();
658     }
659     }
660     } else {
661     // center of mass of the group is the same as position of the atom
662     // if cutoff group does not exist
663     cgConfig->position = config->position;
664 gezelter 1723 cgConfig->velocity = config->velocity;
665 gezelter 1581 }
666    
667 gezelter 1575 fDecomp_->zeroWorkArrays();
668 gezelter 1549 fDecomp_->distributeData();
669 gezelter 1579
670     int cg1, cg2, atom1, atom2, topoDist;
671 gezelter 1723 Vector3d d_grp, dag, d, gvel2, vel2;
672 gezelter 1579 RealType rgrpsq, rgrp, r2, r;
673     RealType electroMult, vdwMult;
674 gezelter 1549 RealType vij;
675 gezelter 1581 Vector3d fij, fg, f1;
676 gezelter 1576 tuple3<RealType, RealType, RealType> cuts;
677 gezelter 1545 RealType rCutSq;
678     bool in_switching_region;
679     RealType sw, dswdr, swderiv;
680 gezelter 1549 vector<int> atomListColumn, atomListRow, atomListLocal;
681 gezelter 1545 InteractionData idat;
682 gezelter 1546 SelfData sdat;
683     RealType mf;
684 gezelter 1579 RealType vpair;
685 jmichalk 1733 RealType dVdFQ1(0.0);
686     RealType dVdFQ2(0.0);
687 gezelter 1583 potVec longRangePotential(0.0);
688     potVec workPot(0.0);
689 gezelter 1760 potVec exPot(0.0);
690 gezelter 1787 Vector3d eField1(0.0);
691     Vector3d eField2(0.0);
692 gezelter 1715 vector<int>::iterator ia, jb;
693 gezelter 1544
694 gezelter 1545 int loopStart, loopEnd;
695 gezelter 1544
696 gezelter 1581 idat.vdwMult = &vdwMult;
697     idat.electroMult = &electroMult;
698 gezelter 1583 idat.pot = &workPot;
699 gezelter 1760 idat.excludedPot = &exPot;
700 gezelter 1583 sdat.pot = fDecomp_->getEmbeddingPotential();
701 gezelter 1761 sdat.excludedPot = fDecomp_->getExcludedSelfPotential();
702 gezelter 1581 idat.vpair = &vpair;
703 jmichalk 1733 idat.dVdFQ1 = &dVdFQ1;
704     idat.dVdFQ2 = &dVdFQ2;
705 gezelter 1787 idat.eField1 = &eField1;
706     idat.eField2 = &eField2;
707 gezelter 1581 idat.f1 = &f1;
708     idat.sw = &sw;
709 gezelter 1583 idat.shiftedPot = (cutoffMethod_ == SHIFTED_POTENTIAL) ? true : false;
710     idat.shiftedForce = (cutoffMethod_ == SHIFTED_FORCE) ? true : false;
711 gezelter 1711 idat.doParticlePot = doParticlePot_;
712 gezelter 1787 idat.doElectricField = doElectricField_;
713 gezelter 1711 sdat.doParticlePot = doParticlePot_;
714 gezelter 1583
715 gezelter 1545 loopEnd = PAIR_LOOP;
716 gezelter 1546 if (info_->requiresPrepair() ) {
717 gezelter 1545 loopStart = PREPAIR_LOOP;
718     } else {
719     loopStart = PAIR_LOOP;
720     }
721 gezelter 1579 for (int iLoop = loopStart; iLoop <= loopEnd; iLoop++) {
722    
723 gezelter 1545 if (iLoop == loopStart) {
724 gezelter 1549 bool update_nlist = fDecomp_->checkNeighborList();
725 gezelter 1545 if (update_nlist)
726 gezelter 1549 neighborList = fDecomp_->buildNeighborList();
727 gezelter 1612 }
728    
729 gezelter 1545 for (vector<pair<int, int> >::iterator it = neighborList.begin();
730     it != neighborList.end(); ++it) {
731 gezelter 1579
732 gezelter 1545 cg1 = (*it).first;
733     cg2 = (*it).second;
734 gezelter 1576
735     cuts = fDecomp_->getGroupCutoffs(cg1, cg2);
736 gezelter 1545
737 gezelter 1549 d_grp = fDecomp_->getIntergroupVector(cg1, cg2);
738 gezelter 1613
739 gezelter 1545 curSnapshot->wrapVector(d_grp);
740     rgrpsq = d_grp.lengthSquare();
741 gezelter 1576 rCutSq = cuts.second;
742    
743 gezelter 1545 if (rgrpsq < rCutSq) {
744 gezelter 1579 idat.rcut = &cuts.first;
745 gezelter 1545 if (iLoop == PAIR_LOOP) {
746 gezelter 1587 vij = 0.0;
747 gezelter 1545 fij = V3Zero;
748 gezelter 1847 eField1 = V3Zero;
749     eField2 = V3Zero;
750 gezelter 1545 }
751    
752 gezelter 1579 in_switching_region = switcher_->getSwitch(rgrpsq, sw, dswdr,
753 gezelter 1576 rgrp);
754 gezelter 1756
755 gezelter 1549 atomListRow = fDecomp_->getAtomsInGroupRow(cg1);
756     atomListColumn = fDecomp_->getAtomsInGroupColumn(cg2);
757 gezelter 1545
758 gezelter 1723 if (doHeatFlux_)
759     gvel2 = fDecomp_->getGroupVelocityColumn(cg2);
760 gezelter 1749
761 gezelter 1715 for (ia = atomListRow.begin();
762 gezelter 1549 ia != atomListRow.end(); ++ia) {
763 gezelter 1545 atom1 = (*ia);
764 gezelter 1749
765 gezelter 1715 for (jb = atomListColumn.begin();
766 gezelter 1549 jb != atomListColumn.end(); ++jb) {
767 gezelter 1545 atom2 = (*jb);
768 gezelter 1593
769 gezelter 1756 if (!fDecomp_->skipAtomPair(atom1, atom2, cg1, cg2)) {
770    
771 gezelter 1579 vpair = 0.0;
772 gezelter 1583 workPot = 0.0;
773 gezelter 1760 exPot = 0.0;
774 gezelter 1581 f1 = V3Zero;
775 jmichalk 1733 dVdFQ1 = 0.0;
776     dVdFQ2 = 0.0;
777 gezelter 1575
778 gezelter 1581 fDecomp_->fillInteractionData(idat, atom1, atom2);
779 gezelter 1749
780 gezelter 1579 topoDist = fDecomp_->getTopologicalDistance(atom1, atom2);
781     vdwMult = vdwScale_[topoDist];
782     electroMult = electrostaticScale_[topoDist];
783 gezelter 1546
784 gezelter 1549 if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
785 gezelter 1579 idat.d = &d_grp;
786     idat.r2 = &rgrpsq;
787 gezelter 1723 if (doHeatFlux_)
788     vel2 = gvel2;
789 gezelter 1545 } else {
790 gezelter 1579 d = fDecomp_->getInteratomicVector(atom1, atom2);
791     curSnapshot->wrapVector( d );
792     r2 = d.lengthSquare();
793     idat.d = &d;
794     idat.r2 = &r2;
795 gezelter 1723 if (doHeatFlux_)
796     vel2 = fDecomp_->getAtomVelocityColumn(atom2);
797 gezelter 1545 }
798 gezelter 1601
799 gezelter 1581 r = sqrt( *(idat.r2) );
800 gezelter 1579 idat.rij = &r;
801 gezelter 1546
802 gezelter 1545 if (iLoop == PREPAIR_LOOP) {
803     interactionMan_->doPrePair(idat);
804     } else {
805     interactionMan_->doPair(idat);
806 gezelter 1575 fDecomp_->unpackInteractionData(idat, atom1, atom2);
807 gezelter 1581 vij += vpair;
808     fij += f1;
809 gezelter 1723 stressTensor -= outProduct( *(idat.d), f1);
810     if (doHeatFlux_)
811     fDecomp_->addToHeatFlux(*(idat.d) * dot(f1, vel2));
812 gezelter 1545 }
813     }
814     }
815     }
816    
817     if (iLoop == PAIR_LOOP) {
818     if (in_switching_region) {
819     swderiv = vij * dswdr / rgrp;
820     fg = swderiv * d_grp;
821     fij += fg;
822    
823 gezelter 1549 if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
824 gezelter 1812 if (!fDecomp_->skipAtomPair(atomListRow[0],
825     atomListColumn[0],
826     cg1, cg2)) {
827     stressTensor -= outProduct( *(idat.d), fg);
828     if (doHeatFlux_)
829     fDecomp_->addToHeatFlux(*(idat.d) * dot(fg, vel2));
830     }
831 gezelter 1545 }
832    
833 gezelter 1715 for (ia = atomListRow.begin();
834 gezelter 1549 ia != atomListRow.end(); ++ia) {
835 gezelter 1545 atom1 = (*ia);
836 gezelter 1569 mf = fDecomp_->getMassFactorRow(atom1);
837 gezelter 1545 // fg is the force on atom ia due to cutoff group's
838     // presence in switching region
839     fg = swderiv * d_grp * mf;
840 gezelter 1549 fDecomp_->addForceToAtomRow(atom1, fg);
841     if (atomListRow.size() > 1) {
842 gezelter 1546 if (info_->usesAtomicVirial()) {
843 gezelter 1545 // find the distance between the atom
844     // and the center of the cutoff group:
845 gezelter 1549 dag = fDecomp_->getAtomToGroupVectorRow(atom1, cg1);
846 gezelter 1723 stressTensor -= outProduct(dag, fg);
847     if (doHeatFlux_)
848     fDecomp_->addToHeatFlux( dag * dot(fg, vel2));
849 gezelter 1545 }
850     }
851     }
852 gezelter 1715 for (jb = atomListColumn.begin();
853 gezelter 1549 jb != atomListColumn.end(); ++jb) {
854 gezelter 1545 atom2 = (*jb);
855 gezelter 1569 mf = fDecomp_->getMassFactorColumn(atom2);
856 gezelter 1545 // fg is the force on atom jb due to cutoff group's
857     // presence in switching region
858     fg = -swderiv * d_grp * mf;
859 gezelter 1549 fDecomp_->addForceToAtomColumn(atom2, fg);
860 gezelter 1545
861 gezelter 1549 if (atomListColumn.size() > 1) {
862 gezelter 1546 if (info_->usesAtomicVirial()) {
863 gezelter 1545 // find the distance between the atom
864     // and the center of the cutoff group:
865 gezelter 1549 dag = fDecomp_->getAtomToGroupVectorColumn(atom2, cg2);
866 gezelter 1723 stressTensor -= outProduct(dag, fg);
867     if (doHeatFlux_)
868     fDecomp_->addToHeatFlux( dag * dot(fg, vel2));
869 gezelter 1545 }
870     }
871     }
872     }
873 gezelter 1613 //if (!info_->usesAtomicVirial()) {
874 gezelter 1723 // stressTensor -= outProduct(d_grp, fij);
875     // if (doHeatFlux_)
876     // fDecomp_->addToHeatFlux( d_grp * dot(fij, vel2));
877 gezelter 1545 //}
878     }
879     }
880     }
881    
882     if (iLoop == PREPAIR_LOOP) {
883 gezelter 1590 if (info_->requiresPrepair()) {
884    
885 gezelter 1549 fDecomp_->collectIntermediateData();
886 gezelter 1570
887 gezelter 1767 for (unsigned int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
888 gezelter 1581 fDecomp_->fillSelfData(sdat, atom1);
889 gezelter 1545 interactionMan_->doPreForce(sdat);
890     }
891 gezelter 1590
892     fDecomp_->distributeIntermediateData();
893    
894 gezelter 1545 }
895     }
896 gezelter 1544 }
897 gezelter 1545
898 gezelter 1756 // collects pairwise information
899 gezelter 1549 fDecomp_->collectData();
900 gezelter 1570
901     if (info_->requiresSelfCorrection()) {
902 gezelter 1767 for (unsigned int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
903 gezelter 1581 fDecomp_->fillSelfData(sdat, atom1);
904 gezelter 1570 interactionMan_->doSelfCorrection(sdat);
905     }
906     }
907    
908 gezelter 1756 // collects single-atom information
909     fDecomp_->collectSelfData();
910    
911 gezelter 1583 longRangePotential = *(fDecomp_->getEmbeddingPotential()) +
912     *(fDecomp_->getPairwisePotential());
913    
914 gezelter 1764 curSnapshot->setLongRangePotential(longRangePotential);
915 gezelter 1794
916     // collects single-atom information
917     fDecomp_->collectSelfData();
918    
919     longRangePotential = *(fDecomp_->getEmbeddingPotential()) +
920     *(fDecomp_->getPairwisePotential());
921    
922     curSnapshot->setLongRangePotential(longRangePotential);
923 gezelter 1761
924     curSnapshot->setExcludedPotentials(*(fDecomp_->getExcludedSelfPotential()) +
925     *(fDecomp_->getExcludedPotential()));
926 gezelter 1760
927 gezelter 507 }
928 gezelter 246
929 gezelter 1126
930 gezelter 1464 void ForceManager::postCalculation() {
931 jmarr 1780
932     vector<Perturbation*>::iterator pi;
933     for (pi = perturbations_.begin(); pi != perturbations_.end(); ++pi) {
934     (*pi)->applyPerturbation();
935     }
936    
937 gezelter 246 SimInfo::MoleculeIterator mi;
938     Molecule* mol;
939     Molecule::RigidBodyIterator rbIter;
940     RigidBody* rb;
941 gezelter 1126 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
942 gezelter 1764
943 gezelter 246 // collect the atomic forces onto rigid bodies
944 gezelter 1126
945     for (mol = info_->beginMolecule(mi); mol != NULL;
946     mol = info_->nextMolecule(mi)) {
947     for (rb = mol->beginRigidBody(rbIter); rb != NULL;
948     rb = mol->nextRigidBody(rbIter)) {
949 gezelter 1464 Mat3x3d rbTau = rb->calcForcesAndTorquesAndVirial();
950 gezelter 1723 stressTensor += rbTau;
951 gezelter 507 }
952 gezelter 1126 }
953 gezelter 1464
954 gezelter 1126 #ifdef IS_MPI
955 gezelter 1723 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, stressTensor.getArrayPointer(), 9,
956     MPI::REALTYPE, MPI::SUM);
957 gezelter 1126 #endif
958 gezelter 1723 curSnapshot->setStressTensor(stressTensor);
959    
960 gezelter 507 }
961 gezelter 1390 } //end namespace OpenMD

Properties

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