ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/ForceManager.cpp
Revision: 2033
Committed: Sat Nov 1 14:12:16 2014 UTC (10 years, 6 months ago) by gezelter
File size: 37656 byte(s)
Log Message:
Fixed a selection issue in ParticleTimeCorrFunc, other whitespace stuff

File Contents

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

Properties

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