ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/ForceManager.cpp
Revision: 1877
Committed: Thu Jun 6 15:43:35 2013 UTC (11 years, 10 months ago) by gezelter
File size: 36930 byte(s)
Log Message:
New electrostatic method, starting to do some performance tuning.

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

Properties

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