ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/ForceManager.cpp
Revision: 1921
Committed: Thu Aug 1 18:23:07 2013 UTC (11 years, 9 months ago) by gezelter
File size: 37288 byte(s)
Log Message:
Fixed a few DumpWriter / FlucQ bugs.  Still working on Ewald, updated 
Quadrupolar test structures.

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/ElectricField.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
398
399
400 void ForceManager::initialize() {
401
402 if (!info_->isTopologyDone()) {
403
404 info_->update();
405 interactionMan_->setSimInfo(info_);
406 interactionMan_->initialize();
407
408 // We want to delay the cutoffs until after the interaction
409 // manager has set up the atom-atom interactions so that we can
410 // query them for suggested cutoff values
411 setupCutoffs();
412
413 info_->prepareTopology();
414
415 doParticlePot_ = info_->getSimParams()->getOutputParticlePotential();
416 doHeatFlux_ = info_->getSimParams()->getPrintHeatFlux();
417 if (doHeatFlux_) doParticlePot_ = true;
418
419 doElectricField_ = info_->getSimParams()->getOutputElectricField();
420
421 }
422
423 ForceFieldOptions& fopts = forceField_->getForceFieldOptions();
424
425 // Force fields can set options on how to scale van der Waals and
426 // electrostatic interactions for atoms connected via bonds, bends
427 // and torsions in this case the topological distance between
428 // atoms is:
429 // 0 = topologically unconnected
430 // 1 = bonded together
431 // 2 = connected via a bend
432 // 3 = connected via a torsion
433
434 vdwScale_.reserve(4);
435 fill(vdwScale_.begin(), vdwScale_.end(), 0.0);
436
437 electrostaticScale_.reserve(4);
438 fill(electrostaticScale_.begin(), electrostaticScale_.end(), 0.0);
439
440 vdwScale_[0] = 1.0;
441 vdwScale_[1] = fopts.getvdw12scale();
442 vdwScale_[2] = fopts.getvdw13scale();
443 vdwScale_[3] = fopts.getvdw14scale();
444
445 electrostaticScale_[0] = 1.0;
446 electrostaticScale_[1] = fopts.getelectrostatic12scale();
447 electrostaticScale_[2] = fopts.getelectrostatic13scale();
448 electrostaticScale_[3] = fopts.getelectrostatic14scale();
449
450 if (info_->getSimParams()->haveElectricField()) {
451 ElectricField* eField = new ElectricField(info_);
452 perturbations_.push_back(eField);
453 }
454
455 usePeriodicBoundaryConditions_ = info_->getSimParams()->getUsePeriodicBoundaryConditions();
456
457 fDecomp_->distributeInitialData();
458
459 initialized_ = true;
460
461 }
462
463 void ForceManager::calcForces() {
464
465 if (!initialized_) initialize();
466
467 preCalculation();
468 shortRangeInteractions();
469 longRangeInteractions();
470 postCalculation();
471 }
472
473 void ForceManager::preCalculation() {
474 SimInfo::MoleculeIterator mi;
475 Molecule* mol;
476 Molecule::AtomIterator ai;
477 Atom* atom;
478 Molecule::RigidBodyIterator rbIter;
479 RigidBody* rb;
480 Molecule::CutoffGroupIterator ci;
481 CutoffGroup* cg;
482
483 // forces and potentials are zeroed here, before any are
484 // accumulated.
485
486 Snapshot* snap = info_->getSnapshotManager()->getCurrentSnapshot();
487
488 snap->setBondPotential(0.0);
489 snap->setBendPotential(0.0);
490 snap->setTorsionPotential(0.0);
491 snap->setInversionPotential(0.0);
492
493 potVec zeroPot(0.0);
494 snap->setLongRangePotential(zeroPot);
495 snap->setExcludedPotentials(zeroPot);
496
497 snap->setRestraintPotential(0.0);
498 snap->setRawPotential(0.0);
499
500 for (mol = info_->beginMolecule(mi); mol != NULL;
501 mol = info_->nextMolecule(mi)) {
502 for(atom = mol->beginAtom(ai); atom != NULL;
503 atom = mol->nextAtom(ai)) {
504 atom->zeroForcesAndTorques();
505 }
506
507 //change the positions of atoms which belong to the rigidbodies
508 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
509 rb = mol->nextRigidBody(rbIter)) {
510 rb->zeroForcesAndTorques();
511 }
512
513 if(info_->getNGlobalCutoffGroups() != info_->getNGlobalAtoms()){
514 for(cg = mol->beginCutoffGroup(ci); cg != NULL;
515 cg = mol->nextCutoffGroup(ci)) {
516 //calculate the center of mass of cutoff group
517 cg->updateCOM();
518 }
519 }
520 }
521
522 // Zero out the stress tensor
523 stressTensor *= 0.0;
524 // Zero out the heatFlux
525 fDecomp_->setHeatFlux( Vector3d(0.0) );
526 }
527
528 void ForceManager::shortRangeInteractions() {
529 Molecule* mol;
530 RigidBody* rb;
531 Bond* bond;
532 Bend* bend;
533 Torsion* torsion;
534 Inversion* inversion;
535 SimInfo::MoleculeIterator mi;
536 Molecule::RigidBodyIterator rbIter;
537 Molecule::BondIterator bondIter;;
538 Molecule::BendIterator bendIter;
539 Molecule::TorsionIterator torsionIter;
540 Molecule::InversionIterator inversionIter;
541 RealType bondPotential = 0.0;
542 RealType bendPotential = 0.0;
543 RealType torsionPotential = 0.0;
544 RealType inversionPotential = 0.0;
545
546 //calculate short range interactions
547 for (mol = info_->beginMolecule(mi); mol != NULL;
548 mol = info_->nextMolecule(mi)) {
549
550 //change the positions of atoms which belong to the rigidbodies
551 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
552 rb = mol->nextRigidBody(rbIter)) {
553 rb->updateAtoms();
554 }
555
556 for (bond = mol->beginBond(bondIter); bond != NULL;
557 bond = mol->nextBond(bondIter)) {
558 bond->calcForce(doParticlePot_);
559 bondPotential += bond->getPotential();
560 }
561
562 for (bend = mol->beginBend(bendIter); bend != NULL;
563 bend = mol->nextBend(bendIter)) {
564
565 RealType angle;
566 bend->calcForce(angle, doParticlePot_);
567 RealType currBendPot = bend->getPotential();
568
569 bendPotential += bend->getPotential();
570 map<Bend*, BendDataSet>::iterator i = bendDataSets.find(bend);
571 if (i == bendDataSets.end()) {
572 BendDataSet dataSet;
573 dataSet.prev.angle = dataSet.curr.angle = angle;
574 dataSet.prev.potential = dataSet.curr.potential = currBendPot;
575 dataSet.deltaV = 0.0;
576 bendDataSets.insert(map<Bend*, BendDataSet>::value_type(bend,
577 dataSet));
578 }else {
579 i->second.prev.angle = i->second.curr.angle;
580 i->second.prev.potential = i->second.curr.potential;
581 i->second.curr.angle = angle;
582 i->second.curr.potential = currBendPot;
583 i->second.deltaV = fabs(i->second.curr.potential -
584 i->second.prev.potential);
585 }
586 }
587
588 for (torsion = mol->beginTorsion(torsionIter); torsion != NULL;
589 torsion = mol->nextTorsion(torsionIter)) {
590 RealType angle;
591 torsion->calcForce(angle, doParticlePot_);
592 RealType currTorsionPot = torsion->getPotential();
593 torsionPotential += torsion->getPotential();
594 map<Torsion*, TorsionDataSet>::iterator i = torsionDataSets.find(torsion);
595 if (i == torsionDataSets.end()) {
596 TorsionDataSet dataSet;
597 dataSet.prev.angle = dataSet.curr.angle = angle;
598 dataSet.prev.potential = dataSet.curr.potential = currTorsionPot;
599 dataSet.deltaV = 0.0;
600 torsionDataSets.insert(map<Torsion*, TorsionDataSet>::value_type(torsion, dataSet));
601 }else {
602 i->second.prev.angle = i->second.curr.angle;
603 i->second.prev.potential = i->second.curr.potential;
604 i->second.curr.angle = angle;
605 i->second.curr.potential = currTorsionPot;
606 i->second.deltaV = fabs(i->second.curr.potential -
607 i->second.prev.potential);
608 }
609 }
610
611 for (inversion = mol->beginInversion(inversionIter);
612 inversion != NULL;
613 inversion = mol->nextInversion(inversionIter)) {
614 RealType angle;
615 inversion->calcForce(angle, doParticlePot_);
616 RealType currInversionPot = inversion->getPotential();
617 inversionPotential += inversion->getPotential();
618 map<Inversion*, InversionDataSet>::iterator i = inversionDataSets.find(inversion);
619 if (i == inversionDataSets.end()) {
620 InversionDataSet dataSet;
621 dataSet.prev.angle = dataSet.curr.angle = angle;
622 dataSet.prev.potential = dataSet.curr.potential = currInversionPot;
623 dataSet.deltaV = 0.0;
624 inversionDataSets.insert(map<Inversion*, InversionDataSet>::value_type(inversion, dataSet));
625 }else {
626 i->second.prev.angle = i->second.curr.angle;
627 i->second.prev.potential = i->second.curr.potential;
628 i->second.curr.angle = angle;
629 i->second.curr.potential = currInversionPot;
630 i->second.deltaV = fabs(i->second.curr.potential -
631 i->second.prev.potential);
632 }
633 }
634 }
635
636 #ifdef IS_MPI
637 // Collect from all nodes. This should eventually be moved into a
638 // SystemDecomposition, but this is a better place than in
639 // Thermo to do the collection.
640 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &bondPotential, 1, MPI::REALTYPE,
641 MPI::SUM);
642 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &bendPotential, 1, MPI::REALTYPE,
643 MPI::SUM);
644 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &torsionPotential, 1,
645 MPI::REALTYPE, MPI::SUM);
646 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &inversionPotential, 1,
647 MPI::REALTYPE, MPI::SUM);
648 #endif
649
650 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
651
652 curSnapshot->setBondPotential(bondPotential);
653 curSnapshot->setBendPotential(bendPotential);
654 curSnapshot->setTorsionPotential(torsionPotential);
655 curSnapshot->setInversionPotential(inversionPotential);
656
657 // RealType shortRangePotential = bondPotential + bendPotential +
658 // torsionPotential + inversionPotential;
659
660 // curSnapshot->setShortRangePotential(shortRangePotential);
661 }
662
663 void ForceManager::longRangeInteractions() {
664
665 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
666 DataStorage* config = &(curSnapshot->atomData);
667 DataStorage* cgConfig = &(curSnapshot->cgData);
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 potVec workPot(0.0);
713 potVec exPot(0.0);
714 Vector3d eField1(0.0);
715 Vector3d eField2(0.0);
716 vector<int>::iterator ia, jb;
717
718 int loopStart, loopEnd;
719
720 idat.rcut = &rCut;
721 idat.vdwMult = &vdwMult;
722 idat.electroMult = &electroMult;
723 idat.pot = &workPot;
724 idat.excludedPot = &exPot;
725 sdat.pot = fDecomp_->getEmbeddingPotential();
726 sdat.excludedPot = fDecomp_->getExcludedSelfPotential();
727 idat.vpair = &vpair;
728 idat.dVdFQ1 = &dVdFQ1;
729 idat.dVdFQ2 = &dVdFQ2;
730 idat.eField1 = &eField1;
731 idat.eField2 = &eField2;
732 idat.f1 = &f1;
733 idat.sw = &sw;
734 idat.shiftedPot = (cutoffMethod_ == SHIFTED_POTENTIAL) ? true : false;
735 idat.shiftedForce = (cutoffMethod_ == SHIFTED_FORCE || cutoffMethod_ == TAYLOR_SHIFTED) ? true : false;
736 idat.doParticlePot = doParticlePot_;
737 idat.doElectricField = doElectricField_;
738 sdat.doParticlePot = doParticlePot_;
739
740 loopEnd = PAIR_LOOP;
741 if (info_->requiresPrepair() ) {
742 loopStart = PREPAIR_LOOP;
743 } else {
744 loopStart = PAIR_LOOP;
745 }
746 for (int iLoop = loopStart; iLoop <= loopEnd; iLoop++) {
747
748 if (iLoop == loopStart) {
749 bool update_nlist = fDecomp_->checkNeighborList();
750 if (update_nlist) {
751 if (!usePeriodicBoundaryConditions_)
752 Mat3x3d bbox = thermo->getBoundingBox();
753 fDecomp_->buildNeighborList(neighborList_);
754 }
755 }
756
757 for (vector<pair<int, int> >::iterator it = neighborList_.begin();
758 it != neighborList_.end(); ++it) {
759
760 cg1 = (*it).first;
761 cg2 = (*it).second;
762
763 fDecomp_->getGroupCutoffs(cg1, cg2, rCut, rCutSq, rListSq);
764
765 d_grp = fDecomp_->getIntergroupVector(cg1, cg2);
766
767 // already wrapped in the getIntergroupVector call:
768 // curSnapshot->wrapVector(d_grp);
769 rgrpsq = d_grp.lengthSquare();
770
771 if (rgrpsq < rCutSq) {
772
773 if (iLoop == PAIR_LOOP) {
774 vij = 0.0;
775 fij.zero();
776 eField1.zero();
777 eField2.zero();
778 }
779
780 in_switching_region = switcher_->getSwitch(rgrpsq, sw, dswdr,
781 rgrp);
782
783 atomListRow = fDecomp_->getAtomsInGroupRow(cg1);
784 atomListColumn = fDecomp_->getAtomsInGroupColumn(cg2);
785
786 if (doHeatFlux_)
787 gvel2 = fDecomp_->getGroupVelocityColumn(cg2);
788
789 for (ia = atomListRow.begin();
790 ia != atomListRow.end(); ++ia) {
791 atom1 = (*ia);
792
793 for (jb = atomListColumn.begin();
794 jb != atomListColumn.end(); ++jb) {
795 atom2 = (*jb);
796
797 if (!fDecomp_->skipAtomPair(atom1, atom2, cg1, cg2)) {
798
799 vpair = 0.0;
800 workPot = 0.0;
801 exPot = 0.0;
802 f1.zero();
803 dVdFQ1 = 0.0;
804 dVdFQ2 = 0.0;
805
806 fDecomp_->fillInteractionData(idat, atom1, atom2);
807
808 topoDist = fDecomp_->getTopologicalDistance(atom1, atom2);
809 vdwMult = vdwScale_[topoDist];
810 electroMult = electrostaticScale_[topoDist];
811
812 if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
813 idat.d = &d_grp;
814 idat.r2 = &rgrpsq;
815 if (doHeatFlux_)
816 vel2 = gvel2;
817 } else {
818 d = fDecomp_->getInteratomicVector(atom1, atom2);
819 curSnapshot->wrapVector( d );
820 r2 = d.lengthSquare();
821 idat.d = &d;
822 idat.r2 = &r2;
823 if (doHeatFlux_)
824 vel2 = fDecomp_->getAtomVelocityColumn(atom2);
825 }
826
827 r = sqrt( *(idat.r2) );
828 idat.rij = &r;
829
830 if (iLoop == PREPAIR_LOOP) {
831 interactionMan_->doPrePair(idat);
832 } else {
833 interactionMan_->doPair(idat);
834 fDecomp_->unpackInteractionData(idat, atom1, atom2);
835 vij += vpair;
836 fij += f1;
837 stressTensor -= outProduct( *(idat.d), f1);
838 if (doHeatFlux_)
839 fDecomp_->addToHeatFlux(*(idat.d) * dot(f1, vel2));
840 }
841 }
842 }
843 }
844
845 if (iLoop == PAIR_LOOP) {
846 if (in_switching_region) {
847 swderiv = vij * dswdr / rgrp;
848 fg = swderiv * d_grp;
849 fij += fg;
850
851 if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
852 if (!fDecomp_->skipAtomPair(atomListRow[0],
853 atomListColumn[0],
854 cg1, cg2)) {
855 stressTensor -= outProduct( *(idat.d), fg);
856 if (doHeatFlux_)
857 fDecomp_->addToHeatFlux(*(idat.d) * dot(fg, vel2));
858 }
859 }
860
861 for (ia = atomListRow.begin();
862 ia != atomListRow.end(); ++ia) {
863 atom1 = (*ia);
864 mf = fDecomp_->getMassFactorRow(atom1);
865 // fg is the force on atom ia due to cutoff group's
866 // presence in switching region
867 fg = swderiv * d_grp * mf;
868 fDecomp_->addForceToAtomRow(atom1, fg);
869 if (atomListRow.size() > 1) {
870 if (info_->usesAtomicVirial()) {
871 // find the distance between the atom
872 // and the center of the cutoff group:
873 dag = fDecomp_->getAtomToGroupVectorRow(atom1, cg1);
874 stressTensor -= outProduct(dag, fg);
875 if (doHeatFlux_)
876 fDecomp_->addToHeatFlux( dag * dot(fg, vel2));
877 }
878 }
879 }
880 for (jb = atomListColumn.begin();
881 jb != atomListColumn.end(); ++jb) {
882 atom2 = (*jb);
883 mf = fDecomp_->getMassFactorColumn(atom2);
884 // fg is the force on atom jb due to cutoff group's
885 // presence in switching region
886 fg = -swderiv * d_grp * mf;
887 fDecomp_->addForceToAtomColumn(atom2, fg);
888
889 if (atomListColumn.size() > 1) {
890 if (info_->usesAtomicVirial()) {
891 // find the distance between the atom
892 // and the center of the cutoff group:
893 dag = fDecomp_->getAtomToGroupVectorColumn(atom2, cg2);
894 stressTensor -= outProduct(dag, fg);
895 if (doHeatFlux_)
896 fDecomp_->addToHeatFlux( dag * dot(fg, vel2));
897 }
898 }
899 }
900 }
901 //if (!info_->usesAtomicVirial()) {
902 // stressTensor -= outProduct(d_grp, fij);
903 // if (doHeatFlux_)
904 // fDecomp_->addToHeatFlux( d_grp * dot(fij, vel2));
905 //}
906 }
907 }
908 }
909
910 if (iLoop == PREPAIR_LOOP) {
911 if (info_->requiresPrepair()) {
912
913 fDecomp_->collectIntermediateData();
914
915 for (unsigned int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
916 fDecomp_->fillSelfData(sdat, atom1);
917 interactionMan_->doPreForce(sdat);
918 }
919
920 fDecomp_->distributeIntermediateData();
921
922 }
923 }
924 }
925
926 // collects pairwise information
927 fDecomp_->collectData();
928 if (cutoffMethod_ == EWALD_FULL) {
929 interactionMan_->doReciprocalSpaceSum(longRangePotential);
930 }
931
932 if (info_->requiresSelfCorrection()) {
933 for (unsigned int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
934 fDecomp_->fillSelfData(sdat, atom1);
935 interactionMan_->doSelfCorrection(sdat);
936 }
937 }
938
939 // collects single-atom information
940 fDecomp_->collectSelfData();
941
942 longRangePotential = *(fDecomp_->getEmbeddingPotential()) +
943 *(fDecomp_->getPairwisePotential());
944
945 curSnapshot->setLongRangePotential(longRangePotential);
946
947 curSnapshot->setExcludedPotentials(*(fDecomp_->getExcludedSelfPotential()) +
948 *(fDecomp_->getExcludedPotential()));
949
950 }
951
952 void ForceManager::postCalculation() {
953
954 vector<Perturbation*>::iterator pi;
955 for (pi = perturbations_.begin(); pi != perturbations_.end(); ++pi) {
956 (*pi)->applyPerturbation();
957 }
958
959 SimInfo::MoleculeIterator mi;
960 Molecule* mol;
961 Molecule::RigidBodyIterator rbIter;
962 RigidBody* rb;
963 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
964
965 // collect the atomic forces onto rigid bodies
966
967 for (mol = info_->beginMolecule(mi); mol != NULL;
968 mol = info_->nextMolecule(mi)) {
969 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
970 rb = mol->nextRigidBody(rbIter)) {
971 Mat3x3d rbTau = rb->calcForcesAndTorquesAndVirial();
972 stressTensor += rbTau;
973 }
974 }
975
976 #ifdef IS_MPI
977 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, stressTensor.getArrayPointer(), 9,
978 MPI::REALTYPE, MPI::SUM);
979 #endif
980 curSnapshot->setStressTensor(stressTensor);
981
982 if (info_->getSimParams()->getUseLongRangeCorrections()) {
983 /*
984 RealType vol = curSnapshot->getVolume();
985 RealType Elrc(0.0);
986 RealType Wlrc(0.0);
987
988 set<AtomType*>::iterator i;
989 set<AtomType*>::iterator j;
990
991 RealType n_i, n_j;
992 RealType rho_i, rho_j;
993 pair<RealType, RealType> LRI;
994
995 for (i = atomTypes_.begin(); i != atomTypes_.end(); ++i) {
996 n_i = RealType(info_->getGlobalCountOfType(*i));
997 rho_i = n_i / vol;
998 for (j = atomTypes_.begin(); j != atomTypes_.end(); ++j) {
999 n_j = RealType(info_->getGlobalCountOfType(*j));
1000 rho_j = n_j / vol;
1001
1002 LRI = interactionMan_->getLongRangeIntegrals( (*i), (*j) );
1003
1004 Elrc += n_i * rho_j * LRI.first;
1005 Wlrc -= rho_i * rho_j * LRI.second;
1006 }
1007 }
1008 Elrc *= 2.0 * NumericConstant::PI;
1009 Wlrc *= 2.0 * NumericConstant::PI;
1010
1011 RealType lrp = curSnapshot->getLongRangePotential();
1012 curSnapshot->setLongRangePotential(lrp + Elrc);
1013 stressTensor += Wlrc * SquareMatrix3<RealType>::identity();
1014 curSnapshot->setStressTensor(stressTensor);
1015 */
1016
1017 }
1018 }
1019 }

Properties

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