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

Properties

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