ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/ForceManager.cpp
Revision: 1780
Committed: Mon Aug 20 18:28:22 2012 UTC (12 years, 8 months ago) by jmarr
File size: 34489 byte(s)
Log Message:
Adding an electric field and the architecture for external perturbations.   Fixing a bug in MultipoleAtomTypesSectionParser.

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 jmarr 1780 #include "perturbations/ElectricField.hpp"
62 gezelter 1549 #include "parallel/ForceMatrixDecomposition.hpp"
63 gezelter 1467
64 gezelter 1583 #include <cstdio>
65     #include <iostream>
66     #include <iomanip>
67    
68 gezelter 1545 using namespace std;
69 gezelter 1390 namespace OpenMD {
70 gezelter 1469
71 gezelter 1545 ForceManager::ForceManager(SimInfo * info) : info_(info) {
72 gezelter 1576 forceField_ = info_->getForceField();
73 gezelter 1577 interactionMan_ = new InteractionManager();
74 gezelter 1579 fDecomp_ = new ForceMatrixDecomposition(info_, interactionMan_);
75 gezelter 1469 }
76 gezelter 1576
77     /**
78     * setupCutoffs
79     *
80 gezelter 1587 * Sets the values of cutoffRadius, switchingRadius, cutoffMethod,
81     * and cutoffPolicy
82 gezelter 1576 *
83     * cutoffRadius : realType
84     * If the cutoffRadius was explicitly set, use that value.
85     * If the cutoffRadius was not explicitly set:
86     * Are there electrostatic atoms? Use 12.0 Angstroms.
87     * No electrostatic atoms? Poll the atom types present in the
88     * simulation for suggested cutoff values (e.g. 2.5 * sigma).
89     * Use the maximum suggested value that was found.
90     *
91 gezelter 1590 * cutoffMethod : (one of HARD, SWITCHED, SHIFTED_FORCE,
92     * or SHIFTED_POTENTIAL)
93 gezelter 1576 * If cutoffMethod was explicitly set, use that choice.
94     * If cutoffMethod was not explicitly set, use SHIFTED_FORCE
95     *
96     * cutoffPolicy : (one of MIX, MAX, TRADITIONAL)
97     * If cutoffPolicy was explicitly set, use that choice.
98     * If cutoffPolicy was not explicitly set, use TRADITIONAL
99 gezelter 1587 *
100     * switchingRadius : realType
101     * If the cutoffMethod was set to SWITCHED:
102     * If the switchingRadius was explicitly set, use that value
103     * (but do a sanity check first).
104     * If the switchingRadius was not explicitly set: use 0.85 *
105     * cutoffRadius_
106     * If the cutoffMethod was not set to SWITCHED:
107     * Set switchingRadius equal to cutoffRadius for safety.
108 gezelter 1576 */
109     void ForceManager::setupCutoffs() {
110    
111     Globals* simParams_ = info_->getSimParams();
112     ForceFieldOptions& forceFieldOptions_ = forceField_->getForceFieldOptions();
113 gezelter 1613 int mdFileVersion;
114 jmichalk 1754 rCut_ = 0.0; //Needs a value for a later max() call;
115 gezelter 1755
116 gezelter 1613 if (simParams_->haveMDfileVersion())
117     mdFileVersion = simParams_->getMDfileVersion();
118     else
119     mdFileVersion = 0;
120    
121 gezelter 1576 if (simParams_->haveCutoffRadius()) {
122     rCut_ = simParams_->getCutoffRadius();
123     } else {
124     if (info_->usesElectrostaticAtoms()) {
125     sprintf(painCave.errMsg,
126     "ForceManager::setupCutoffs: No value was set for the cutoffRadius.\n"
127     "\tOpenMD will use a default value of 12.0 angstroms"
128     "\tfor the cutoffRadius.\n");
129     painCave.isFatal = 0;
130     painCave.severity = OPENMD_INFO;
131     simError();
132     rCut_ = 12.0;
133     } else {
134     RealType thisCut;
135     set<AtomType*>::iterator i;
136     set<AtomType*> atomTypes;
137     atomTypes = info_->getSimulatedAtomTypes();
138     for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
139     thisCut = interactionMan_->getSuggestedCutoffRadius((*i));
140     rCut_ = max(thisCut, rCut_);
141     }
142     sprintf(painCave.errMsg,
143     "ForceManager::setupCutoffs: No value was set for the cutoffRadius.\n"
144     "\tOpenMD will use %lf angstroms.\n",
145     rCut_);
146     painCave.isFatal = 0;
147     painCave.severity = OPENMD_INFO;
148     simError();
149 gezelter 1579 }
150 gezelter 1576 }
151    
152 gezelter 1583 fDecomp_->setUserCutoff(rCut_);
153 gezelter 1584 interactionMan_->setCutoffRadius(rCut_);
154 gezelter 1583
155 gezelter 1576 map<string, CutoffMethod> stringToCutoffMethod;
156     stringToCutoffMethod["HARD"] = HARD;
157     stringToCutoffMethod["SWITCHED"] = SWITCHED;
158     stringToCutoffMethod["SHIFTED_POTENTIAL"] = SHIFTED_POTENTIAL;
159     stringToCutoffMethod["SHIFTED_FORCE"] = SHIFTED_FORCE;
160 gezelter 1545
161 gezelter 1576 if (simParams_->haveCutoffMethod()) {
162     string cutMeth = toUpperCopy(simParams_->getCutoffMethod());
163     map<string, CutoffMethod>::iterator i;
164     i = stringToCutoffMethod.find(cutMeth);
165     if (i == stringToCutoffMethod.end()) {
166     sprintf(painCave.errMsg,
167     "ForceManager::setupCutoffs: Could not find chosen cutoffMethod %s\n"
168     "\tShould be one of: "
169     "HARD, SWITCHED, SHIFTED_POTENTIAL, or SHIFTED_FORCE\n",
170     cutMeth.c_str());
171     painCave.isFatal = 1;
172     painCave.severity = OPENMD_ERROR;
173     simError();
174     } else {
175     cutoffMethod_ = i->second;
176     }
177     } else {
178 gezelter 1616 if (mdFileVersion > 1) {
179     sprintf(painCave.errMsg,
180     "ForceManager::setupCutoffs: No value was set for the cutoffMethod.\n"
181     "\tOpenMD will use SHIFTED_FORCE.\n");
182     painCave.isFatal = 0;
183     painCave.severity = OPENMD_INFO;
184     simError();
185     cutoffMethod_ = SHIFTED_FORCE;
186     } else {
187     // handle the case where the old file version was in play
188     // (there should be no cutoffMethod, so we have to deduce it
189     // from other data).
190    
191     sprintf(painCave.errMsg,
192     "ForceManager::setupCutoffs : DEPRECATED FILE FORMAT!\n"
193     "\tOpenMD found a file which does not set a cutoffMethod.\n"
194     "\tOpenMD will attempt to deduce a cutoffMethod using the\n"
195     "\tbehavior of the older (version 1) code. To remove this\n"
196     "\twarning, add an explicit cutoffMethod and change the top\n"
197     "\tof the file so that it begins with <OpenMD version=2>\n");
198     painCave.isFatal = 0;
199     painCave.severity = OPENMD_WARNING;
200     simError();
201    
202     // The old file version tethered the shifting behavior to the
203     // electrostaticSummationMethod keyword.
204    
205     if (simParams_->haveElectrostaticSummationMethod()) {
206 gezelter 1710 string myMethod = simParams_->getElectrostaticSummationMethod();
207 gezelter 1616 toUpper(myMethod);
208    
209     if (myMethod == "SHIFTED_POTENTIAL") {
210     cutoffMethod_ = SHIFTED_POTENTIAL;
211     } else if (myMethod == "SHIFTED_FORCE") {
212     cutoffMethod_ = SHIFTED_FORCE;
213     }
214    
215     if (simParams_->haveSwitchingRadius())
216     rSwitch_ = simParams_->getSwitchingRadius();
217    
218     if (myMethod == "SHIFTED_POTENTIAL" || myMethod == "SHIFTED_FORCE") {
219     if (simParams_->haveSwitchingRadius()){
220     sprintf(painCave.errMsg,
221     "ForceManager::setupCutoffs : DEPRECATED ERROR MESSAGE\n"
222     "\tA value was set for the switchingRadius\n"
223     "\teven though the electrostaticSummationMethod was\n"
224     "\tset to %s\n", myMethod.c_str());
225     painCave.severity = OPENMD_WARNING;
226     painCave.isFatal = 1;
227     simError();
228     }
229     }
230     if (abs(rCut_ - rSwitch_) < 0.0001) {
231     if (cutoffMethod_ == SHIFTED_FORCE) {
232     sprintf(painCave.errMsg,
233     "ForceManager::setupCutoffs : DEPRECATED BEHAVIOR\n"
234     "\tcutoffRadius and switchingRadius are set to the\n"
235     "\tsame value. OpenMD will use shifted force\n"
236     "\tpotentials instead of switching functions.\n");
237     painCave.isFatal = 0;
238     painCave.severity = OPENMD_WARNING;
239     simError();
240     } else {
241     cutoffMethod_ = SHIFTED_POTENTIAL;
242     sprintf(painCave.errMsg,
243     "ForceManager::setupCutoffs : DEPRECATED BEHAVIOR\n"
244     "\tcutoffRadius and switchingRadius are set to the\n"
245     "\tsame value. OpenMD will use shifted potentials\n"
246     "\tinstead of switching functions.\n");
247     painCave.isFatal = 0;
248     painCave.severity = OPENMD_WARNING;
249     simError();
250     }
251     }
252     }
253     }
254 gezelter 1576 }
255    
256     map<string, CutoffPolicy> stringToCutoffPolicy;
257     stringToCutoffPolicy["MIX"] = MIX;
258     stringToCutoffPolicy["MAX"] = MAX;
259     stringToCutoffPolicy["TRADITIONAL"] = TRADITIONAL;
260    
261 gezelter 1710 string cutPolicy;
262 gezelter 1576 if (forceFieldOptions_.haveCutoffPolicy()){
263     cutPolicy = forceFieldOptions_.getCutoffPolicy();
264     }else if (simParams_->haveCutoffPolicy()) {
265     cutPolicy = simParams_->getCutoffPolicy();
266     }
267    
268     if (!cutPolicy.empty()){
269     toUpper(cutPolicy);
270     map<string, CutoffPolicy>::iterator i;
271     i = stringToCutoffPolicy.find(cutPolicy);
272    
273     if (i == stringToCutoffPolicy.end()) {
274     sprintf(painCave.errMsg,
275     "ForceManager::setupCutoffs: Could not find chosen cutoffPolicy %s\n"
276     "\tShould be one of: "
277     "MIX, MAX, or TRADITIONAL\n",
278     cutPolicy.c_str());
279     painCave.isFatal = 1;
280     painCave.severity = OPENMD_ERROR;
281     simError();
282     } else {
283     cutoffPolicy_ = i->second;
284     }
285     } else {
286     sprintf(painCave.errMsg,
287     "ForceManager::setupCutoffs: No value was set for the cutoffPolicy.\n"
288     "\tOpenMD will use TRADITIONAL.\n");
289     painCave.isFatal = 0;
290     painCave.severity = OPENMD_INFO;
291     simError();
292     cutoffPolicy_ = TRADITIONAL;
293     }
294 gezelter 1587
295 gezelter 1579 fDecomp_->setCutoffPolicy(cutoffPolicy_);
296 gezelter 1587
297     // create the switching function object:
298 gezelter 1576
299 gezelter 1577 switcher_ = new SwitchingFunction();
300 gezelter 1587
301     if (cutoffMethod_ == SWITCHED) {
302     if (simParams_->haveSwitchingRadius()) {
303     rSwitch_ = simParams_->getSwitchingRadius();
304     if (rSwitch_ > rCut_) {
305     sprintf(painCave.errMsg,
306     "ForceManager::setupCutoffs: switchingRadius (%f) is larger "
307     "than the cutoffRadius(%f)\n", rSwitch_, rCut_);
308     painCave.isFatal = 1;
309     painCave.severity = OPENMD_ERROR;
310     simError();
311     }
312     } else {
313     rSwitch_ = 0.85 * rCut_;
314 gezelter 1576 sprintf(painCave.errMsg,
315 gezelter 1587 "ForceManager::setupCutoffs: No value was set for the switchingRadius.\n"
316     "\tOpenMD will use a default value of 85 percent of the cutoffRadius.\n"
317     "\tswitchingRadius = %f. for this simulation\n", rSwitch_);
318     painCave.isFatal = 0;
319     painCave.severity = OPENMD_WARNING;
320 gezelter 1576 simError();
321     }
322 gezelter 1587 } else {
323 gezelter 1618 if (mdFileVersion > 1) {
324     // throw an error if we define a switching radius and don't need one.
325     // older file versions should not do this.
326     if (simParams_->haveSwitchingRadius()) {
327     map<string, CutoffMethod>::const_iterator it;
328     string theMeth;
329     for (it = stringToCutoffMethod.begin();
330     it != stringToCutoffMethod.end(); ++it) {
331     if (it->second == cutoffMethod_) {
332     theMeth = it->first;
333     break;
334     }
335 gezelter 1587 }
336 gezelter 1618 sprintf(painCave.errMsg,
337     "ForceManager::setupCutoffs: the cutoffMethod (%s)\n"
338     "\tis not set to SWITCHED, so switchingRadius value\n"
339     "\twill be ignored for this simulation\n", theMeth.c_str());
340     painCave.isFatal = 0;
341     painCave.severity = OPENMD_WARNING;
342     simError();
343 gezelter 1587 }
344     }
345     rSwitch_ = rCut_;
346     }
347 gezelter 1576
348 gezelter 1577 // Default to cubic switching function.
349     sft_ = cubic;
350 gezelter 1576 if (simParams_->haveSwitchingFunctionType()) {
351     string funcType = simParams_->getSwitchingFunctionType();
352     toUpper(funcType);
353     if (funcType == "CUBIC") {
354     sft_ = cubic;
355     } else {
356     if (funcType == "FIFTH_ORDER_POLYNOMIAL") {
357     sft_ = fifth_order_poly;
358     } else {
359     // throw error
360     sprintf( painCave.errMsg,
361     "ForceManager::setupSwitching : Unknown switchingFunctionType. (Input file specified %s .)\n"
362     "\tswitchingFunctionType must be one of: "
363     "\"cubic\" or \"fifth_order_polynomial\".",
364     funcType.c_str() );
365     painCave.isFatal = 1;
366     painCave.severity = OPENMD_ERROR;
367     simError();
368     }
369     }
370     }
371     switcher_->setSwitchType(sft_);
372     switcher_->setSwitch(rSwitch_, rCut_);
373 gezelter 1584 interactionMan_->setSwitchingRadius(rSwitch_);
374 gezelter 1576 }
375 gezelter 1616
376    
377    
378 gezelter 1576
379     void ForceManager::initialize() {
380    
381 gezelter 1569 if (!info_->isTopologyDone()) {
382 gezelter 1590
383 gezelter 507 info_->update();
384 gezelter 1546 interactionMan_->setSimInfo(info_);
385     interactionMan_->initialize();
386 gezelter 1576
387     // We want to delay the cutoffs until after the interaction
388     // manager has set up the atom-atom interactions so that we can
389     // query them for suggested cutoff values
390     setupCutoffs();
391    
392     info_->prepareTopology();
393 gezelter 1711
394     doParticlePot_ = info_->getSimParams()->getOutputParticlePotential();
395 gezelter 1723 doHeatFlux_ = info_->getSimParams()->getPrintHeatFlux();
396     if (doHeatFlux_) doParticlePot_ = true;
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 1715 vector<int>::iterator ia, jb;
691 gezelter 1544
692 gezelter 1545 int loopStart, loopEnd;
693 gezelter 1544
694 gezelter 1581 idat.vdwMult = &vdwMult;
695     idat.electroMult = &electroMult;
696 gezelter 1583 idat.pot = &workPot;
697 gezelter 1760 idat.excludedPot = &exPot;
698 gezelter 1583 sdat.pot = fDecomp_->getEmbeddingPotential();
699 gezelter 1761 sdat.excludedPot = fDecomp_->getExcludedSelfPotential();
700 gezelter 1581 idat.vpair = &vpair;
701 jmichalk 1733 idat.dVdFQ1 = &dVdFQ1;
702     idat.dVdFQ2 = &dVdFQ2;
703 gezelter 1581 idat.f1 = &f1;
704     idat.sw = &sw;
705 gezelter 1583 idat.shiftedPot = (cutoffMethod_ == SHIFTED_POTENTIAL) ? true : false;
706     idat.shiftedForce = (cutoffMethod_ == SHIFTED_FORCE) ? true : false;
707 gezelter 1711 idat.doParticlePot = doParticlePot_;
708     sdat.doParticlePot = doParticlePot_;
709 gezelter 1583
710 gezelter 1545 loopEnd = PAIR_LOOP;
711 gezelter 1546 if (info_->requiresPrepair() ) {
712 gezelter 1545 loopStart = PREPAIR_LOOP;
713     } else {
714     loopStart = PAIR_LOOP;
715     }
716 gezelter 1579 for (int iLoop = loopStart; iLoop <= loopEnd; iLoop++) {
717    
718 gezelter 1545 if (iLoop == loopStart) {
719 gezelter 1549 bool update_nlist = fDecomp_->checkNeighborList();
720 gezelter 1545 if (update_nlist)
721 gezelter 1549 neighborList = fDecomp_->buildNeighborList();
722 gezelter 1612 }
723    
724 gezelter 1545 for (vector<pair<int, int> >::iterator it = neighborList.begin();
725     it != neighborList.end(); ++it) {
726 gezelter 1579
727 gezelter 1545 cg1 = (*it).first;
728     cg2 = (*it).second;
729 gezelter 1576
730     cuts = fDecomp_->getGroupCutoffs(cg1, cg2);
731 gezelter 1545
732 gezelter 1549 d_grp = fDecomp_->getIntergroupVector(cg1, cg2);
733 gezelter 1613
734 gezelter 1545 curSnapshot->wrapVector(d_grp);
735     rgrpsq = d_grp.lengthSquare();
736 gezelter 1576 rCutSq = cuts.second;
737    
738 gezelter 1545 if (rgrpsq < rCutSq) {
739 gezelter 1579 idat.rcut = &cuts.first;
740 gezelter 1545 if (iLoop == PAIR_LOOP) {
741 gezelter 1587 vij = 0.0;
742 gezelter 1545 fij = V3Zero;
743     }
744    
745 gezelter 1579 in_switching_region = switcher_->getSwitch(rgrpsq, sw, dswdr,
746 gezelter 1576 rgrp);
747 gezelter 1756
748 gezelter 1549 atomListRow = fDecomp_->getAtomsInGroupRow(cg1);
749     atomListColumn = fDecomp_->getAtomsInGroupColumn(cg2);
750 gezelter 1545
751 gezelter 1723 if (doHeatFlux_)
752     gvel2 = fDecomp_->getGroupVelocityColumn(cg2);
753 gezelter 1749
754 gezelter 1715 for (ia = atomListRow.begin();
755 gezelter 1549 ia != atomListRow.end(); ++ia) {
756 gezelter 1545 atom1 = (*ia);
757 gezelter 1749
758 gezelter 1715 for (jb = atomListColumn.begin();
759 gezelter 1549 jb != atomListColumn.end(); ++jb) {
760 gezelter 1545 atom2 = (*jb);
761 gezelter 1593
762 gezelter 1756 if (!fDecomp_->skipAtomPair(atom1, atom2, cg1, cg2)) {
763    
764 gezelter 1579 vpair = 0.0;
765 gezelter 1583 workPot = 0.0;
766 gezelter 1760 exPot = 0.0;
767 gezelter 1581 f1 = V3Zero;
768 jmichalk 1733 dVdFQ1 = 0.0;
769     dVdFQ2 = 0.0;
770 gezelter 1575
771 gezelter 1581 fDecomp_->fillInteractionData(idat, atom1, atom2);
772 gezelter 1749
773 gezelter 1579 topoDist = fDecomp_->getTopologicalDistance(atom1, atom2);
774     vdwMult = vdwScale_[topoDist];
775     electroMult = electrostaticScale_[topoDist];
776 gezelter 1546
777 gezelter 1549 if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
778 gezelter 1579 idat.d = &d_grp;
779     idat.r2 = &rgrpsq;
780 gezelter 1723 if (doHeatFlux_)
781     vel2 = gvel2;
782 gezelter 1545 } else {
783 gezelter 1579 d = fDecomp_->getInteratomicVector(atom1, atom2);
784     curSnapshot->wrapVector( d );
785     r2 = d.lengthSquare();
786     idat.d = &d;
787     idat.r2 = &r2;
788 gezelter 1723 if (doHeatFlux_)
789     vel2 = fDecomp_->getAtomVelocityColumn(atom2);
790 gezelter 1545 }
791 gezelter 1601
792 gezelter 1581 r = sqrt( *(idat.r2) );
793 gezelter 1579 idat.rij = &r;
794 gezelter 1546
795 gezelter 1545 if (iLoop == PREPAIR_LOOP) {
796     interactionMan_->doPrePair(idat);
797     } else {
798     interactionMan_->doPair(idat);
799 gezelter 1575 fDecomp_->unpackInteractionData(idat, atom1, atom2);
800 gezelter 1581 vij += vpair;
801     fij += f1;
802 gezelter 1723 stressTensor -= outProduct( *(idat.d), f1);
803     if (doHeatFlux_)
804     fDecomp_->addToHeatFlux(*(idat.d) * dot(f1, vel2));
805 gezelter 1545 }
806     }
807     }
808     }
809    
810     if (iLoop == PAIR_LOOP) {
811     if (in_switching_region) {
812     swderiv = vij * dswdr / rgrp;
813     fg = swderiv * d_grp;
814     fij += fg;
815    
816 gezelter 1549 if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
817 gezelter 1723 stressTensor -= outProduct( *(idat.d), fg);
818     if (doHeatFlux_)
819     fDecomp_->addToHeatFlux(*(idat.d) * dot(fg, vel2));
820    
821 gezelter 1545 }
822    
823 gezelter 1715 for (ia = atomListRow.begin();
824 gezelter 1549 ia != atomListRow.end(); ++ia) {
825 gezelter 1545 atom1 = (*ia);
826 gezelter 1569 mf = fDecomp_->getMassFactorRow(atom1);
827 gezelter 1545 // fg is the force on atom ia due to cutoff group's
828     // presence in switching region
829     fg = swderiv * d_grp * mf;
830 gezelter 1549 fDecomp_->addForceToAtomRow(atom1, fg);
831     if (atomListRow.size() > 1) {
832 gezelter 1546 if (info_->usesAtomicVirial()) {
833 gezelter 1545 // find the distance between the atom
834     // and the center of the cutoff group:
835 gezelter 1549 dag = fDecomp_->getAtomToGroupVectorRow(atom1, cg1);
836 gezelter 1723 stressTensor -= outProduct(dag, fg);
837     if (doHeatFlux_)
838     fDecomp_->addToHeatFlux( dag * dot(fg, vel2));
839 gezelter 1545 }
840     }
841     }
842 gezelter 1715 for (jb = atomListColumn.begin();
843 gezelter 1549 jb != atomListColumn.end(); ++jb) {
844 gezelter 1545 atom2 = (*jb);
845 gezelter 1569 mf = fDecomp_->getMassFactorColumn(atom2);
846 gezelter 1545 // fg is the force on atom jb due to cutoff group's
847     // presence in switching region
848     fg = -swderiv * d_grp * mf;
849 gezelter 1549 fDecomp_->addForceToAtomColumn(atom2, fg);
850 gezelter 1545
851 gezelter 1549 if (atomListColumn.size() > 1) {
852 gezelter 1546 if (info_->usesAtomicVirial()) {
853 gezelter 1545 // find the distance between the atom
854     // and the center of the cutoff group:
855 gezelter 1549 dag = fDecomp_->getAtomToGroupVectorColumn(atom2, cg2);
856 gezelter 1723 stressTensor -= outProduct(dag, fg);
857     if (doHeatFlux_)
858     fDecomp_->addToHeatFlux( dag * dot(fg, vel2));
859 gezelter 1545 }
860     }
861     }
862     }
863 gezelter 1613 //if (!info_->usesAtomicVirial()) {
864 gezelter 1723 // stressTensor -= outProduct(d_grp, fij);
865     // if (doHeatFlux_)
866     // fDecomp_->addToHeatFlux( d_grp * dot(fij, vel2));
867 gezelter 1545 //}
868     }
869     }
870     }
871    
872     if (iLoop == PREPAIR_LOOP) {
873 gezelter 1590 if (info_->requiresPrepair()) {
874    
875 gezelter 1549 fDecomp_->collectIntermediateData();
876 gezelter 1570
877 gezelter 1767 for (unsigned int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
878 gezelter 1581 fDecomp_->fillSelfData(sdat, atom1);
879 gezelter 1545 interactionMan_->doPreForce(sdat);
880     }
881 gezelter 1590
882     fDecomp_->distributeIntermediateData();
883    
884 gezelter 1545 }
885     }
886 gezelter 1544 }
887 gezelter 1545
888 gezelter 1756 // collects pairwise information
889 gezelter 1549 fDecomp_->collectData();
890 gezelter 1570
891     if (info_->requiresSelfCorrection()) {
892 gezelter 1767 for (unsigned int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
893 gezelter 1581 fDecomp_->fillSelfData(sdat, atom1);
894 gezelter 1570 interactionMan_->doSelfCorrection(sdat);
895     }
896     }
897    
898 gezelter 1756 // collects single-atom information
899     fDecomp_->collectSelfData();
900    
901 gezelter 1583 longRangePotential = *(fDecomp_->getEmbeddingPotential()) +
902     *(fDecomp_->getPairwisePotential());
903    
904 gezelter 1764 curSnapshot->setLongRangePotential(longRangePotential);
905 gezelter 1761
906     curSnapshot->setExcludedPotentials(*(fDecomp_->getExcludedSelfPotential()) +
907     *(fDecomp_->getExcludedPotential()));
908 gezelter 1760
909 gezelter 507 }
910 gezelter 246
911 gezelter 1126
912 gezelter 1464 void ForceManager::postCalculation() {
913 jmarr 1780
914     vector<Perturbation*>::iterator pi;
915     for (pi = perturbations_.begin(); pi != perturbations_.end(); ++pi) {
916     (*pi)->applyPerturbation();
917     }
918    
919 gezelter 246 SimInfo::MoleculeIterator mi;
920     Molecule* mol;
921     Molecule::RigidBodyIterator rbIter;
922     RigidBody* rb;
923 gezelter 1126 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
924 gezelter 1764
925 gezelter 246 // collect the atomic forces onto rigid bodies
926 gezelter 1126
927     for (mol = info_->beginMolecule(mi); mol != NULL;
928     mol = info_->nextMolecule(mi)) {
929     for (rb = mol->beginRigidBody(rbIter); rb != NULL;
930     rb = mol->nextRigidBody(rbIter)) {
931 gezelter 1464 Mat3x3d rbTau = rb->calcForcesAndTorquesAndVirial();
932 gezelter 1723 stressTensor += rbTau;
933 gezelter 507 }
934 gezelter 1126 }
935 gezelter 1464
936 gezelter 1126 #ifdef IS_MPI
937 gezelter 1723 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, stressTensor.getArrayPointer(), 9,
938     MPI::REALTYPE, MPI::SUM);
939 gezelter 1126 #endif
940 gezelter 1723 curSnapshot->setStressTensor(stressTensor);
941    
942 gezelter 507 }
943 gezelter 1390 } //end namespace OpenMD

Properties

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