ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/ForceManager.cpp
Revision: 1896
Committed: Tue Jul 2 20:02:31 2013 UTC (11 years, 10 months ago) by gezelter
File size: 36987 byte(s)
Log Message:
Speedup tests

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     #include "perturbations/ElectricField.hpp"
61     #include "parallel/ForceMatrixDecomposition.hpp"
62    
63     #include <cstdio>
64     #include <iostream>
65     #include <iomanip>
66    
67     using namespace std;
68 gezelter 1390 namespace OpenMD {
69 gezelter 1782
70 gezelter 1879 ForceManager::ForceManager(SimInfo * info) : info_(info), switcher_(NULL),
71     initialized_(false) {
72 gezelter 1782 forceField_ = info_->getForceField();
73     interactionMan_ = new InteractionManager();
74     fDecomp_ = new ForceMatrixDecomposition(info_, interactionMan_);
75 gezelter 1879 thermo = new Thermo(info_);
76 gezelter 1782 }
77 gezelter 246
78 gezelter 1879 ForceManager::~ForceManager() {
79     perturbations_.clear();
80    
81     delete switcher_;
82     delete interactionMan_;
83     delete fDecomp_;
84     delete thermo;
85     }
86    
87 gezelter 1782 /**
88     * setupCutoffs
89     *
90     * Sets the values of cutoffRadius, switchingRadius, cutoffMethod,
91     * and cutoffPolicy
92     *
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 1879 * cutoffMethod : (one of HARD, SWITCHED, SHIFTED_FORCE, TAYLOR_SHIFTED,
102 gezelter 1782 * or SHIFTED_POTENTIAL)
103     * 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     *
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     */
119     void ForceManager::setupCutoffs() {
120 gezelter 1126
121 gezelter 1782 Globals* simParams_ = info_->getSimParams();
122     ForceFieldOptions& forceFieldOptions_ = forceField_->getForceFieldOptions();
123     int mdFileVersion;
124     rCut_ = 0.0; //Needs a value for a later max() call;
125    
126     if (simParams_->haveMDfileVersion())
127     mdFileVersion = simParams_->getMDfileVersion();
128     else
129     mdFileVersion = 0;
130    
131 gezelter 1879 // 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 1782 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 1879 for (i = atomTypes_.begin(); i != atomTypes_.end(); ++i) {
153 gezelter 1782 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     }
164     }
165    
166     fDecomp_->setUserCutoff(rCut_);
167     interactionMan_->setCutoffRadius(rCut_);
168    
169     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 1879 stringToCutoffMethod["TAYLOR_SHIFTED"] = TAYLOR_SHIFTED;
175 gezelter 1782
176     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 1879 "HARD, SWITCHED, SHIFTED_POTENTIAL, TAYLOR_SHIFTED,\n"
185     "\tor SHIFTED_FORCE\n",
186 gezelter 1782 cutMeth.c_str());
187     painCave.isFatal = 1;
188     painCave.severity = OPENMD_ERROR;
189     simError();
190     } else {
191     cutoffMethod_ = i->second;
192     }
193     } else {
194     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     string myMethod = simParams_->getElectrostaticSummationMethod();
223     toUpper(myMethod);
224    
225     if (myMethod == "SHIFTED_POTENTIAL") {
226     cutoffMethod_ = SHIFTED_POTENTIAL;
227     } else if (myMethod == "SHIFTED_FORCE") {
228     cutoffMethod_ = SHIFTED_FORCE;
229 gezelter 1879 } else if (myMethod == "TAYLOR_SHIFTED") {
230     cutoffMethod_ = TAYLOR_SHIFTED;
231 gezelter 1782 }
232    
233     if (simParams_->haveSwitchingRadius())
234     rSwitch_ = simParams_->getSwitchingRadius();
235    
236 gezelter 1879 if (myMethod == "SHIFTED_POTENTIAL" || myMethod == "SHIFTED_FORCE" ||
237     myMethod == "TAYLOR_SHIFTED") {
238 gezelter 1782 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     }
274    
275     map<string, CutoffPolicy> stringToCutoffPolicy;
276     stringToCutoffPolicy["MIX"] = MIX;
277     stringToCutoffPolicy["MAX"] = MAX;
278     stringToCutoffPolicy["TRADITIONAL"] = TRADITIONAL;
279    
280     string cutPolicy;
281     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    
314     fDecomp_->setCutoffPolicy(cutoffPolicy_);
315    
316     // create the switching function object:
317    
318     switcher_ = new SwitchingFunction();
319    
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     sprintf(painCave.errMsg,
334     "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     simError();
340     }
341     } else {
342     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     }
355     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     }
363     }
364     rSwitch_ = rCut_;
365     }
366    
367     // Default to cubic switching function.
368     sft_ = cubic;
369     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    
394    
395    
396    
397     void ForceManager::initialize() {
398    
399     if (!info_->isTopologyDone()) {
400    
401 gezelter 507 info_->update();
402 gezelter 1782 interactionMan_->setSimInfo(info_);
403     interactionMan_->initialize();
404    
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    
412     doParticlePot_ = info_->getSimParams()->getOutputParticlePotential();
413     doHeatFlux_ = info_->getSimParams()->getPrintHeatFlux();
414     if (doHeatFlux_) doParticlePot_ = true;
415 gezelter 1879
416     doElectricField_ = info_->getSimParams()->getOutputElectricField();
417 gezelter 1782
418 gezelter 246 }
419 gezelter 1782
420     ForceFieldOptions& fopts = forceField_->getForceFieldOptions();
421 gezelter 1126
422 gezelter 1782 // 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     // 0 = topologically unconnected
427     // 1 = bonded together
428     // 2 = connected via a bend
429     // 3 = connected via a torsion
430 gezelter 246
431 gezelter 1782 vdwScale_.reserve(4);
432     fill(vdwScale_.begin(), vdwScale_.end(), 0.0);
433 gezelter 246
434 gezelter 1782 electrostaticScale_.reserve(4);
435     fill(electrostaticScale_.begin(), electrostaticScale_.end(), 0.0);
436 gezelter 246
437 gezelter 1782 vdwScale_[0] = 1.0;
438     vdwScale_[1] = fopts.getvdw12scale();
439     vdwScale_[2] = fopts.getvdw13scale();
440     vdwScale_[3] = fopts.getvdw14scale();
441 tim 749
442 gezelter 1782 electrostaticScale_[0] = 1.0;
443     electrostaticScale_[1] = fopts.getelectrostatic12scale();
444     electrostaticScale_[2] = fopts.getelectrostatic13scale();
445     electrostaticScale_[3] = fopts.getelectrostatic14scale();
446    
447     if (info_->getSimParams()->haveElectricField()) {
448     ElectricField* eField = new ElectricField(info_);
449     perturbations_.push_back(eField);
450     }
451    
452 gezelter 1879 usePeriodicBoundaryConditions_ = info_->getSimParams()->getUsePeriodicBoundaryConditions();
453    
454 gezelter 1782 fDecomp_->distributeInitialData();
455 gezelter 1879
456 gezelter 1782 initialized_ = true;
457 gezelter 1879
458 gezelter 507 }
459 gezelter 1879
460 gezelter 1782 void ForceManager::calcForces() {
461    
462     if (!initialized_) initialize();
463 gezelter 1879
464 gezelter 1782 preCalculation();
465     shortRangeInteractions();
466     longRangeInteractions();
467     postCalculation();
468     }
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 1782 Molecule::CutoffGroupIterator ci;
478     CutoffGroup* cg;
479 gezelter 246
480 gezelter 1782 // forces and potentials are zeroed here, before any are
481     // accumulated.
482 chuckv 1245
483 gezelter 1782 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 1782 for(atom = mol->beginAtom(ai); atom != NULL;
500     atom = mol->nextAtom(ai)) {
501 gezelter 507 atom->zeroForcesAndTorques();
502     }
503 gezelter 1782
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 1782
510     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    
519 gezelter 1126 // Zero out the stress tensor
520 gezelter 1782 stressTensor *= 0.0;
521     // Zero out the heatFlux
522     fDecomp_->setHeatFlux( Vector3d(0.0) );
523 gezelter 507 }
524 gezelter 1126
525 gezelter 1782 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 1782 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 1782 bend->calcForce(angle, doParticlePot_);
564 gezelter 1126 RealType currBendPot = bend->getPotential();
565 gezelter 1448
566 gezelter 1126 bendPotential += bend->getPotential();
567 gezelter 1782 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 1782 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 1782 torsion->calcForce(angle, doParticlePot_);
589 tim 963 RealType currTorsionPot = torsion->getPotential();
590 gezelter 1126 torsionPotential += torsion->getPotential();
591 gezelter 1782 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 1782 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 1782
608 cli2 1275 for (inversion = mol->beginInversion(inversionIter);
609     inversion != NULL;
610     inversion = mol->nextInversion(inversionIter)) {
611     RealType angle;
612 gezelter 1782 inversion->calcForce(angle, doParticlePot_);
613 cli2 1275 RealType currInversionPot = inversion->getPotential();
614     inversionPotential += inversion->getPotential();
615 gezelter 1782 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 1782 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 1782
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 gezelter 246 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
648 gezelter 1782
649     curSnapshot->setBondPotential(bondPotential);
650     curSnapshot->setBendPotential(bendPotential);
651     curSnapshot->setTorsionPotential(torsionPotential);
652     curSnapshot->setInversionPotential(inversionPotential);
653 tim 665
654 gezelter 1782 // RealType shortRangePotential = bondPotential + bendPotential +
655     // torsionPotential + inversionPotential;
656    
657     // curSnapshot->setShortRangePotential(shortRangePotential);
658 gezelter 507 }
659 gezelter 1126
660 gezelter 1782 void ForceManager::longRangeInteractions() {
661 gezelter 246
662 gezelter 1782 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
663     DataStorage* config = &(curSnapshot->atomData);
664     DataStorage* cgConfig = &(curSnapshot->cgData);
665    
666 gezelter 246 //calculate the center of mass of cutoff group
667 gezelter 1782
668 gezelter 246 SimInfo::MoleculeIterator mi;
669     Molecule* mol;
670     Molecule::CutoffGroupIterator ci;
671     CutoffGroup* cg;
672 gezelter 1782
673     if(info_->getNCutoffGroups() > 0){
674 gezelter 1126 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 gezelter 1782 cg->updateCOM();
679 gezelter 246 }
680 gezelter 1782 }
681 gezelter 246 } else {
682 gezelter 1126 // center of mass of the group is the same as position of the atom
683     // if cutoff group does not exist
684 gezelter 1782 cgConfig->position = config->position;
685     cgConfig->velocity = config->velocity;
686 gezelter 246 }
687 gezelter 1782
688     fDecomp_->zeroWorkArrays();
689     fDecomp_->distributeData();
690 gezelter 1126
691 gezelter 1782 int cg1, cg2, atom1, atom2, topoDist;
692     Vector3d d_grp, dag, d, gvel2, vel2;
693     RealType rgrpsq, rgrp, r2, r;
694     RealType electroMult, vdwMult;
695     RealType vij;
696     Vector3d fij, fg, f1;
697     tuple3<RealType, RealType, RealType> cuts;
698 gezelter 1896 RealType rCut, rCutSq, rListSq;
699 gezelter 1782 bool in_switching_region;
700     RealType sw, dswdr, swderiv;
701 gezelter 1879 vector<int> atomListColumn, atomListRow;
702 gezelter 1782 InteractionData idat;
703     SelfData sdat;
704     RealType mf;
705     RealType vpair;
706     RealType dVdFQ1(0.0);
707     RealType dVdFQ2(0.0);
708     potVec longRangePotential(0.0);
709     potVec workPot(0.0);
710     potVec exPot(0.0);
711 gezelter 1880 Vector3d eField1(0.0);
712     Vector3d eField2(0.0);
713 gezelter 1782 vector<int>::iterator ia, jb;
714 gezelter 246
715 gezelter 1782 int loopStart, loopEnd;
716 gezelter 1895
717 gezelter 1896 idat.rcut = &rCut;
718 gezelter 1782 idat.vdwMult = &vdwMult;
719     idat.electroMult = &electroMult;
720     idat.pot = &workPot;
721     idat.excludedPot = &exPot;
722     sdat.pot = fDecomp_->getEmbeddingPotential();
723     sdat.excludedPot = fDecomp_->getExcludedSelfPotential();
724     idat.vpair = &vpair;
725     idat.dVdFQ1 = &dVdFQ1;
726     idat.dVdFQ2 = &dVdFQ2;
727 gezelter 1788 idat.eField1 = &eField1;
728 gezelter 1879 idat.eField2 = &eField2;
729 gezelter 1782 idat.f1 = &f1;
730     idat.sw = &sw;
731     idat.shiftedPot = (cutoffMethod_ == SHIFTED_POTENTIAL) ? true : false;
732 gezelter 1879 idat.shiftedForce = (cutoffMethod_ == SHIFTED_FORCE || cutoffMethod_ == TAYLOR_SHIFTED) ? true : false;
733 gezelter 1782 idat.doParticlePot = doParticlePot_;
734 gezelter 1879 idat.doElectricField = doElectricField_;
735 gezelter 1782 sdat.doParticlePot = doParticlePot_;
736    
737     loopEnd = PAIR_LOOP;
738     if (info_->requiresPrepair() ) {
739     loopStart = PREPAIR_LOOP;
740     } else {
741     loopStart = PAIR_LOOP;
742 chuckv 664 }
743 gezelter 1782 for (int iLoop = loopStart; iLoop <= loopEnd; iLoop++) {
744 gezelter 1126
745 gezelter 1782 if (iLoop == loopStart) {
746     bool update_nlist = fDecomp_->checkNeighborList();
747 gezelter 1879 if (update_nlist) {
748     if (!usePeriodicBoundaryConditions_)
749     Mat3x3d bbox = thermo->getBoundingBox();
750 gezelter 1895 fDecomp_->buildNeighborList(neighborList_);
751 gezelter 1879 }
752     }
753 gezelter 1782
754 gezelter 1895 for (vector<pair<int, int> >::iterator it = neighborList_.begin();
755     it != neighborList_.end(); ++it) {
756 gezelter 1782
757     cg1 = (*it).first;
758     cg2 = (*it).second;
759    
760 gezelter 1896 fDecomp_->getGroupCutoffs(cg1, cg2, rCut, rCutSq, rListSq);
761 gezelter 1782
762     d_grp = fDecomp_->getIntergroupVector(cg1, cg2);
763    
764 gezelter 1893 // already wrapped in the getIntergroupVector call:
765     // curSnapshot->wrapVector(d_grp);
766 gezelter 1782 rgrpsq = d_grp.lengthSquare();
767    
768     if (rgrpsq < rCutSq) {
769 gezelter 1896
770 gezelter 1782 if (iLoop == PAIR_LOOP) {
771     vij = 0.0;
772 gezelter 1879 fij.zero();
773     eField1.zero();
774     eField2.zero();
775 gezelter 1782 }
776    
777     in_switching_region = switcher_->getSwitch(rgrpsq, sw, dswdr,
778     rgrp);
779    
780     atomListRow = fDecomp_->getAtomsInGroupRow(cg1);
781     atomListColumn = fDecomp_->getAtomsInGroupColumn(cg2);
782    
783     if (doHeatFlux_)
784     gvel2 = fDecomp_->getGroupVelocityColumn(cg2);
785    
786     for (ia = atomListRow.begin();
787     ia != atomListRow.end(); ++ia) {
788     atom1 = (*ia);
789    
790     for (jb = atomListColumn.begin();
791     jb != atomListColumn.end(); ++jb) {
792     atom2 = (*jb);
793    
794     if (!fDecomp_->skipAtomPair(atom1, atom2, cg1, cg2)) {
795    
796     vpair = 0.0;
797     workPot = 0.0;
798     exPot = 0.0;
799 gezelter 1879 f1.zero();
800 gezelter 1782 dVdFQ1 = 0.0;
801     dVdFQ2 = 0.0;
802    
803     fDecomp_->fillInteractionData(idat, atom1, atom2);
804    
805     topoDist = fDecomp_->getTopologicalDistance(atom1, atom2);
806     vdwMult = vdwScale_[topoDist];
807     electroMult = electrostaticScale_[topoDist];
808    
809     if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
810     idat.d = &d_grp;
811     idat.r2 = &rgrpsq;
812     if (doHeatFlux_)
813     vel2 = gvel2;
814     } else {
815     d = fDecomp_->getInteratomicVector(atom1, atom2);
816     curSnapshot->wrapVector( d );
817     r2 = d.lengthSquare();
818     idat.d = &d;
819     idat.r2 = &r2;
820     if (doHeatFlux_)
821     vel2 = fDecomp_->getAtomVelocityColumn(atom2);
822     }
823    
824     r = sqrt( *(idat.r2) );
825     idat.rij = &r;
826    
827     if (iLoop == PREPAIR_LOOP) {
828     interactionMan_->doPrePair(idat);
829     } else {
830     interactionMan_->doPair(idat);
831     fDecomp_->unpackInteractionData(idat, atom1, atom2);
832     vij += vpair;
833     fij += f1;
834     stressTensor -= outProduct( *(idat.d), f1);
835     if (doHeatFlux_)
836     fDecomp_->addToHeatFlux(*(idat.d) * dot(f1, vel2));
837     }
838     }
839     }
840     }
841    
842     if (iLoop == PAIR_LOOP) {
843     if (in_switching_region) {
844     swderiv = vij * dswdr / rgrp;
845     fg = swderiv * d_grp;
846     fij += fg;
847    
848     if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
849 gezelter 1809 if (!fDecomp_->skipAtomPair(atomListRow[0],
850     atomListColumn[0],
851     cg1, cg2)) {
852     stressTensor -= outProduct( *(idat.d), fg);
853     if (doHeatFlux_)
854     fDecomp_->addToHeatFlux(*(idat.d) * dot(fg, vel2));
855     }
856 gezelter 1782 }
857    
858     for (ia = atomListRow.begin();
859     ia != atomListRow.end(); ++ia) {
860     atom1 = (*ia);
861     mf = fDecomp_->getMassFactorRow(atom1);
862     // fg is the force on atom ia due to cutoff group's
863     // presence in switching region
864     fg = swderiv * d_grp * mf;
865     fDecomp_->addForceToAtomRow(atom1, fg);
866     if (atomListRow.size() > 1) {
867     if (info_->usesAtomicVirial()) {
868     // find the distance between the atom
869     // and the center of the cutoff group:
870     dag = fDecomp_->getAtomToGroupVectorRow(atom1, cg1);
871     stressTensor -= outProduct(dag, fg);
872     if (doHeatFlux_)
873     fDecomp_->addToHeatFlux( dag * dot(fg, vel2));
874     }
875     }
876     }
877     for (jb = atomListColumn.begin();
878     jb != atomListColumn.end(); ++jb) {
879     atom2 = (*jb);
880     mf = fDecomp_->getMassFactorColumn(atom2);
881     // fg is the force on atom jb due to cutoff group's
882     // presence in switching region
883     fg = -swderiv * d_grp * mf;
884     fDecomp_->addForceToAtomColumn(atom2, fg);
885    
886     if (atomListColumn.size() > 1) {
887     if (info_->usesAtomicVirial()) {
888     // find the distance between the atom
889     // and the center of the cutoff group:
890     dag = fDecomp_->getAtomToGroupVectorColumn(atom2, cg2);
891     stressTensor -= outProduct(dag, fg);
892     if (doHeatFlux_)
893     fDecomp_->addToHeatFlux( dag * dot(fg, vel2));
894     }
895     }
896     }
897     }
898     //if (!info_->usesAtomicVirial()) {
899     // stressTensor -= outProduct(d_grp, fij);
900     // if (doHeatFlux_)
901     // fDecomp_->addToHeatFlux( d_grp * dot(fij, vel2));
902     //}
903     }
904     }
905     }
906    
907     if (iLoop == PREPAIR_LOOP) {
908     if (info_->requiresPrepair()) {
909    
910     fDecomp_->collectIntermediateData();
911    
912     for (unsigned int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
913     fDecomp_->fillSelfData(sdat, atom1);
914     interactionMan_->doPreForce(sdat);
915     }
916    
917     fDecomp_->distributeIntermediateData();
918    
919     }
920     }
921 gezelter 246 }
922 gezelter 1880
923 gezelter 1782 // collects pairwise information
924     fDecomp_->collectData();
925    
926     if (info_->requiresSelfCorrection()) {
927     for (unsigned int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
928     fDecomp_->fillSelfData(sdat, atom1);
929     interactionMan_->doSelfCorrection(sdat);
930     }
931 chrisfen 998 }
932 gezelter 1782
933     // collects single-atom information
934     fDecomp_->collectSelfData();
935    
936     longRangePotential = *(fDecomp_->getEmbeddingPotential()) +
937     *(fDecomp_->getPairwisePotential());
938    
939     curSnapshot->setLongRangePotential(longRangePotential);
940 gezelter 1126
941 gezelter 1782 curSnapshot->setExcludedPotentials(*(fDecomp_->getExcludedSelfPotential()) +
942     *(fDecomp_->getExcludedPotential()));
943    
944 gezelter 507 }
945 gezelter 246
946 gezelter 1126
947 gezelter 1464 void ForceManager::postCalculation() {
948 gezelter 1782
949     vector<Perturbation*>::iterator pi;
950     for (pi = perturbations_.begin(); pi != perturbations_.end(); ++pi) {
951     (*pi)->applyPerturbation();
952     }
953    
954 gezelter 246 SimInfo::MoleculeIterator mi;
955     Molecule* mol;
956     Molecule::RigidBodyIterator rbIter;
957     RigidBody* rb;
958 gezelter 1126 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
959 gezelter 1782
960 gezelter 246 // collect the atomic forces onto rigid bodies
961 gezelter 1126
962     for (mol = info_->beginMolecule(mi); mol != NULL;
963     mol = info_->nextMolecule(mi)) {
964     for (rb = mol->beginRigidBody(rbIter); rb != NULL;
965     rb = mol->nextRigidBody(rbIter)) {
966 gezelter 1464 Mat3x3d rbTau = rb->calcForcesAndTorquesAndVirial();
967 gezelter 1782 stressTensor += rbTau;
968 gezelter 507 }
969 gezelter 1126 }
970 gezelter 1464
971 gezelter 1126 #ifdef IS_MPI
972 gezelter 1782 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, stressTensor.getArrayPointer(), 9,
973     MPI::REALTYPE, MPI::SUM);
974 gezelter 1126 #endif
975 gezelter 1782 curSnapshot->setStressTensor(stressTensor);
976    
977 gezelter 1879 if (info_->getSimParams()->getUseLongRangeCorrections()) {
978     /*
979     RealType vol = curSnapshot->getVolume();
980     RealType Elrc(0.0);
981     RealType Wlrc(0.0);
982    
983     set<AtomType*>::iterator i;
984     set<AtomType*>::iterator j;
985    
986     RealType n_i, n_j;
987     RealType rho_i, rho_j;
988     pair<RealType, RealType> LRI;
989    
990     for (i = atomTypes_.begin(); i != atomTypes_.end(); ++i) {
991     n_i = RealType(info_->getGlobalCountOfType(*i));
992     rho_i = n_i / vol;
993     for (j = atomTypes_.begin(); j != atomTypes_.end(); ++j) {
994     n_j = RealType(info_->getGlobalCountOfType(*j));
995     rho_j = n_j / vol;
996    
997     LRI = interactionMan_->getLongRangeIntegrals( (*i), (*j) );
998    
999     Elrc += n_i * rho_j * LRI.first;
1000     Wlrc -= rho_i * rho_j * LRI.second;
1001     }
1002     }
1003     Elrc *= 2.0 * NumericConstant::PI;
1004     Wlrc *= 2.0 * NumericConstant::PI;
1005    
1006     RealType lrp = curSnapshot->getLongRangePotential();
1007     curSnapshot->setLongRangePotential(lrp + Elrc);
1008     stressTensor += Wlrc * SquareMatrix3<RealType>::identity();
1009     curSnapshot->setStressTensor(stressTensor);
1010     */
1011    
1012     }
1013 gezelter 507 }
1014 gezelter 1879 }

Properties

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