ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/ForceManager.cpp
Revision: 1868
Committed: Tue Apr 30 15:56:54 2013 UTC (12 years ago) by gezelter
File size: 36593 byte(s)
Log Message:
CLearing out some memory leaks

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

Properties

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