ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/ForceManager.cpp
Revision: 1874
Committed: Wed May 15 15:09:35 2013 UTC (11 years, 11 months ago) by gezelter
File size: 36646 byte(s)
Log Message:
Fixed a bunch of cppcheck warnings.

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

Properties

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