ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/ForceManager.cpp
Revision: 1993
Committed: Tue Apr 29 17:32:31 2014 UTC (11 years ago) by gezelter
File size: 37611 byte(s)
Log Message:
Added sitePotentials for the Choi et al. potential-frequency maps for nitriles

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

Properties

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