ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/ForceManager.cpp
Revision: 1744
Committed: Tue Jun 5 18:07:08 2012 UTC (12 years, 10 months ago) by gezelter
File size: 32982 byte(s)
Log Message:
Fixes for minimization

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, 24107 (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 * @time 10:39am
48 * @version 1.0
49 */
50
51
52 #include "brains/ForceManager.hpp"
53 #include "primitives/Molecule.hpp"
54 #define __OPENMD_C
55 #include "utils/simError.h"
56 #include "primitives/Bond.hpp"
57 #include "primitives/Bend.hpp"
58 #include "primitives/Torsion.hpp"
59 #include "primitives/Inversion.hpp"
60 #include "nonbonded/NonBondedInteraction.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) {
71 forceField_ = info_->getForceField();
72 interactionMan_ = new InteractionManager();
73 fDecomp_ = new ForceMatrixDecomposition(info_, interactionMan_);
74 }
75
76 /**
77 * setupCutoffs
78 *
79 * Sets the values of cutoffRadius, switchingRadius, cutoffMethod,
80 * and cutoffPolicy
81 *
82 * cutoffRadius : realType
83 * If the cutoffRadius was explicitly set, use that value.
84 * If the cutoffRadius was not explicitly set:
85 * Are there electrostatic atoms? Use 12.0 Angstroms.
86 * No electrostatic atoms? Poll the atom types present in the
87 * simulation for suggested cutoff values (e.g. 2.5 * sigma).
88 * Use the maximum suggested value that was found.
89 *
90 * cutoffMethod : (one of HARD, SWITCHED, SHIFTED_FORCE,
91 * or SHIFTED_POTENTIAL)
92 * If cutoffMethod was explicitly set, use that choice.
93 * If cutoffMethod was not explicitly set, use SHIFTED_FORCE
94 *
95 * cutoffPolicy : (one of MIX, MAX, TRADITIONAL)
96 * If cutoffPolicy was explicitly set, use that choice.
97 * If cutoffPolicy was not explicitly set, use TRADITIONAL
98 *
99 * switchingRadius : realType
100 * If the cutoffMethod was set to SWITCHED:
101 * If the switchingRadius was explicitly set, use that value
102 * (but do a sanity check first).
103 * If the switchingRadius was not explicitly set: use 0.85 *
104 * cutoffRadius_
105 * If the cutoffMethod was not set to SWITCHED:
106 * Set switchingRadius equal to cutoffRadius for safety.
107 */
108 void ForceManager::setupCutoffs() {
109
110 Globals* simParams_ = info_->getSimParams();
111 ForceFieldOptions& forceFieldOptions_ = forceField_->getForceFieldOptions();
112 int mdFileVersion;
113
114 if (simParams_->haveMDfileVersion())
115 mdFileVersion = simParams_->getMDfileVersion();
116 else
117 mdFileVersion = 0;
118
119 if (simParams_->haveCutoffRadius()) {
120 rCut_ = simParams_->getCutoffRadius();
121 } else {
122 if (info_->usesElectrostaticAtoms()) {
123 sprintf(painCave.errMsg,
124 "ForceManager::setupCutoffs: No value was set for the cutoffRadius.\n"
125 "\tOpenMD will use a default value of 12.0 angstroms"
126 "\tfor the cutoffRadius.\n");
127 painCave.isFatal = 0;
128 painCave.severity = OPENMD_INFO;
129 simError();
130 rCut_ = 12.0;
131 } else {
132 RealType thisCut;
133 set<AtomType*>::iterator i;
134 set<AtomType*> atomTypes;
135 atomTypes = info_->getSimulatedAtomTypes();
136 for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
137 thisCut = interactionMan_->getSuggestedCutoffRadius((*i));
138 rCut_ = max(thisCut, rCut_);
139 }
140 sprintf(painCave.errMsg,
141 "ForceManager::setupCutoffs: No value was set for the cutoffRadius.\n"
142 "\tOpenMD will use %lf angstroms.\n",
143 rCut_);
144 painCave.isFatal = 0;
145 painCave.severity = OPENMD_INFO;
146 simError();
147 }
148 }
149
150 fDecomp_->setUserCutoff(rCut_);
151 interactionMan_->setCutoffRadius(rCut_);
152
153 map<string, CutoffMethod> stringToCutoffMethod;
154 stringToCutoffMethod["HARD"] = HARD;
155 stringToCutoffMethod["SWITCHED"] = SWITCHED;
156 stringToCutoffMethod["SHIFTED_POTENTIAL"] = SHIFTED_POTENTIAL;
157 stringToCutoffMethod["SHIFTED_FORCE"] = SHIFTED_FORCE;
158
159 if (simParams_->haveCutoffMethod()) {
160 string cutMeth = toUpperCopy(simParams_->getCutoffMethod());
161 map<string, CutoffMethod>::iterator i;
162 i = stringToCutoffMethod.find(cutMeth);
163 if (i == stringToCutoffMethod.end()) {
164 sprintf(painCave.errMsg,
165 "ForceManager::setupCutoffs: Could not find chosen cutoffMethod %s\n"
166 "\tShould be one of: "
167 "HARD, SWITCHED, SHIFTED_POTENTIAL, or SHIFTED_FORCE\n",
168 cutMeth.c_str());
169 painCave.isFatal = 1;
170 painCave.severity = OPENMD_ERROR;
171 simError();
172 } else {
173 cutoffMethod_ = i->second;
174 }
175 } else {
176 if (mdFileVersion > 1) {
177 sprintf(painCave.errMsg,
178 "ForceManager::setupCutoffs: No value was set for the cutoffMethod.\n"
179 "\tOpenMD will use SHIFTED_FORCE.\n");
180 painCave.isFatal = 0;
181 painCave.severity = OPENMD_INFO;
182 simError();
183 cutoffMethod_ = SHIFTED_FORCE;
184 } else {
185 // handle the case where the old file version was in play
186 // (there should be no cutoffMethod, so we have to deduce it
187 // from other data).
188
189 sprintf(painCave.errMsg,
190 "ForceManager::setupCutoffs : DEPRECATED FILE FORMAT!\n"
191 "\tOpenMD found a file which does not set a cutoffMethod.\n"
192 "\tOpenMD will attempt to deduce a cutoffMethod using the\n"
193 "\tbehavior of the older (version 1) code. To remove this\n"
194 "\twarning, add an explicit cutoffMethod and change the top\n"
195 "\tof the file so that it begins with <OpenMD version=2>\n");
196 painCave.isFatal = 0;
197 painCave.severity = OPENMD_WARNING;
198 simError();
199
200 // The old file version tethered the shifting behavior to the
201 // electrostaticSummationMethod keyword.
202
203 if (simParams_->haveElectrostaticSummationMethod()) {
204 string myMethod = simParams_->getElectrostaticSummationMethod();
205 toUpper(myMethod);
206
207 if (myMethod == "SHIFTED_POTENTIAL") {
208 cutoffMethod_ = SHIFTED_POTENTIAL;
209 } else if (myMethod == "SHIFTED_FORCE") {
210 cutoffMethod_ = SHIFTED_FORCE;
211 }
212
213 if (simParams_->haveSwitchingRadius())
214 rSwitch_ = simParams_->getSwitchingRadius();
215
216 if (myMethod == "SHIFTED_POTENTIAL" || myMethod == "SHIFTED_FORCE") {
217 if (simParams_->haveSwitchingRadius()){
218 sprintf(painCave.errMsg,
219 "ForceManager::setupCutoffs : DEPRECATED ERROR MESSAGE\n"
220 "\tA value was set for the switchingRadius\n"
221 "\teven though the electrostaticSummationMethod was\n"
222 "\tset to %s\n", myMethod.c_str());
223 painCave.severity = OPENMD_WARNING;
224 painCave.isFatal = 1;
225 simError();
226 }
227 }
228 if (abs(rCut_ - rSwitch_) < 0.0001) {
229 if (cutoffMethod_ == SHIFTED_FORCE) {
230 sprintf(painCave.errMsg,
231 "ForceManager::setupCutoffs : DEPRECATED BEHAVIOR\n"
232 "\tcutoffRadius and switchingRadius are set to the\n"
233 "\tsame value. OpenMD will use shifted force\n"
234 "\tpotentials instead of switching functions.\n");
235 painCave.isFatal = 0;
236 painCave.severity = OPENMD_WARNING;
237 simError();
238 } else {
239 cutoffMethod_ = SHIFTED_POTENTIAL;
240 sprintf(painCave.errMsg,
241 "ForceManager::setupCutoffs : DEPRECATED BEHAVIOR\n"
242 "\tcutoffRadius and switchingRadius are set to the\n"
243 "\tsame value. OpenMD will use shifted potentials\n"
244 "\tinstead of switching functions.\n");
245 painCave.isFatal = 0;
246 painCave.severity = OPENMD_WARNING;
247 simError();
248 }
249 }
250 }
251 }
252 }
253
254 map<string, CutoffPolicy> stringToCutoffPolicy;
255 stringToCutoffPolicy["MIX"] = MIX;
256 stringToCutoffPolicy["MAX"] = MAX;
257 stringToCutoffPolicy["TRADITIONAL"] = TRADITIONAL;
258
259 string cutPolicy;
260 if (forceFieldOptions_.haveCutoffPolicy()){
261 cutPolicy = forceFieldOptions_.getCutoffPolicy();
262 }else if (simParams_->haveCutoffPolicy()) {
263 cutPolicy = simParams_->getCutoffPolicy();
264 }
265
266 if (!cutPolicy.empty()){
267 toUpper(cutPolicy);
268 map<string, CutoffPolicy>::iterator i;
269 i = stringToCutoffPolicy.find(cutPolicy);
270
271 if (i == stringToCutoffPolicy.end()) {
272 sprintf(painCave.errMsg,
273 "ForceManager::setupCutoffs: Could not find chosen cutoffPolicy %s\n"
274 "\tShould be one of: "
275 "MIX, MAX, or TRADITIONAL\n",
276 cutPolicy.c_str());
277 painCave.isFatal = 1;
278 painCave.severity = OPENMD_ERROR;
279 simError();
280 } else {
281 cutoffPolicy_ = i->second;
282 }
283 } else {
284 sprintf(painCave.errMsg,
285 "ForceManager::setupCutoffs: No value was set for the cutoffPolicy.\n"
286 "\tOpenMD will use TRADITIONAL.\n");
287 painCave.isFatal = 0;
288 painCave.severity = OPENMD_INFO;
289 simError();
290 cutoffPolicy_ = TRADITIONAL;
291 }
292
293 fDecomp_->setCutoffPolicy(cutoffPolicy_);
294
295 // create the switching function object:
296
297 switcher_ = new SwitchingFunction();
298
299 if (cutoffMethod_ == SWITCHED) {
300 if (simParams_->haveSwitchingRadius()) {
301 rSwitch_ = simParams_->getSwitchingRadius();
302 if (rSwitch_ > rCut_) {
303 sprintf(painCave.errMsg,
304 "ForceManager::setupCutoffs: switchingRadius (%f) is larger "
305 "than the cutoffRadius(%f)\n", rSwitch_, rCut_);
306 painCave.isFatal = 1;
307 painCave.severity = OPENMD_ERROR;
308 simError();
309 }
310 } else {
311 rSwitch_ = 0.85 * rCut_;
312 sprintf(painCave.errMsg,
313 "ForceManager::setupCutoffs: No value was set for the switchingRadius.\n"
314 "\tOpenMD will use a default value of 85 percent of the cutoffRadius.\n"
315 "\tswitchingRadius = %f. for this simulation\n", rSwitch_);
316 painCave.isFatal = 0;
317 painCave.severity = OPENMD_WARNING;
318 simError();
319 }
320 } else {
321 if (mdFileVersion > 1) {
322 // throw an error if we define a switching radius and don't need one.
323 // older file versions should not do this.
324 if (simParams_->haveSwitchingRadius()) {
325 map<string, CutoffMethod>::const_iterator it;
326 string theMeth;
327 for (it = stringToCutoffMethod.begin();
328 it != stringToCutoffMethod.end(); ++it) {
329 if (it->second == cutoffMethod_) {
330 theMeth = it->first;
331 break;
332 }
333 }
334 sprintf(painCave.errMsg,
335 "ForceManager::setupCutoffs: the cutoffMethod (%s)\n"
336 "\tis not set to SWITCHED, so switchingRadius value\n"
337 "\twill be ignored for this simulation\n", theMeth.c_str());
338 painCave.isFatal = 0;
339 painCave.severity = OPENMD_WARNING;
340 simError();
341 }
342 }
343 rSwitch_ = rCut_;
344 }
345
346 // Default to cubic switching function.
347 sft_ = cubic;
348 if (simParams_->haveSwitchingFunctionType()) {
349 string funcType = simParams_->getSwitchingFunctionType();
350 toUpper(funcType);
351 if (funcType == "CUBIC") {
352 sft_ = cubic;
353 } else {
354 if (funcType == "FIFTH_ORDER_POLYNOMIAL") {
355 sft_ = fifth_order_poly;
356 } else {
357 // throw error
358 sprintf( painCave.errMsg,
359 "ForceManager::setupSwitching : Unknown switchingFunctionType. (Input file specified %s .)\n"
360 "\tswitchingFunctionType must be one of: "
361 "\"cubic\" or \"fifth_order_polynomial\".",
362 funcType.c_str() );
363 painCave.isFatal = 1;
364 painCave.severity = OPENMD_ERROR;
365 simError();
366 }
367 }
368 }
369 switcher_->setSwitchType(sft_);
370 switcher_->setSwitch(rSwitch_, rCut_);
371 interactionMan_->setSwitchingRadius(rSwitch_);
372 }
373
374
375
376
377 void ForceManager::initialize() {
378
379 if (!info_->isTopologyDone()) {
380
381 info_->update();
382 interactionMan_->setSimInfo(info_);
383 interactionMan_->initialize();
384
385 // We want to delay the cutoffs until after the interaction
386 // manager has set up the atom-atom interactions so that we can
387 // query them for suggested cutoff values
388 setupCutoffs();
389
390 info_->prepareTopology();
391
392 doParticlePot_ = info_->getSimParams()->getOutputParticlePotential();
393 doHeatFlux_ = info_->getSimParams()->getPrintHeatFlux();
394 if (doHeatFlux_) doParticlePot_ = true;
395
396 }
397
398 ForceFieldOptions& fopts = forceField_->getForceFieldOptions();
399
400 // Force fields can set options on how to scale van der Waals and
401 // electrostatic interactions for atoms connected via bonds, bends
402 // and torsions in this case the topological distance between
403 // atoms is:
404 // 0 = topologically unconnected
405 // 1 = bonded together
406 // 2 = connected via a bend
407 // 3 = connected via a torsion
408
409 vdwScale_.reserve(4);
410 fill(vdwScale_.begin(), vdwScale_.end(), 0.0);
411
412 electrostaticScale_.reserve(4);
413 fill(electrostaticScale_.begin(), electrostaticScale_.end(), 0.0);
414
415 vdwScale_[0] = 1.0;
416 vdwScale_[1] = fopts.getvdw12scale();
417 vdwScale_[2] = fopts.getvdw13scale();
418 vdwScale_[3] = fopts.getvdw14scale();
419
420 electrostaticScale_[0] = 1.0;
421 electrostaticScale_[1] = fopts.getelectrostatic12scale();
422 electrostaticScale_[2] = fopts.getelectrostatic13scale();
423 electrostaticScale_[3] = fopts.getelectrostatic14scale();
424
425 fDecomp_->distributeInitialData();
426
427 initialized_ = true;
428
429 }
430
431 void ForceManager::calcForces() {
432
433 if (!initialized_) initialize();
434
435 preCalculation();
436 shortRangeInteractions();
437 longRangeInteractions();
438 postCalculation();
439 }
440
441 void ForceManager::preCalculation() {
442 SimInfo::MoleculeIterator mi;
443 Molecule* mol;
444 Molecule::AtomIterator ai;
445 Atom* atom;
446 Molecule::RigidBodyIterator rbIter;
447 RigidBody* rb;
448 Molecule::CutoffGroupIterator ci;
449 CutoffGroup* cg;
450
451 // forces are zeroed here, before any are accumulated.
452
453 for (mol = info_->beginMolecule(mi); mol != NULL;
454 mol = info_->nextMolecule(mi)) {
455 for(atom = mol->beginAtom(ai); atom != NULL;
456 atom = mol->nextAtom(ai)) {
457 atom->zeroForcesAndTorques();
458 }
459
460 //change the positions of atoms which belong to the rigidbodies
461 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
462 rb = mol->nextRigidBody(rbIter)) {
463 rb->zeroForcesAndTorques();
464 }
465
466 if(info_->getNGlobalCutoffGroups() != info_->getNGlobalAtoms()){
467 for(cg = mol->beginCutoffGroup(ci); cg != NULL;
468 cg = mol->nextCutoffGroup(ci)) {
469 //calculate the center of mass of cutoff group
470 cg->updateCOM();
471 }
472 }
473 }
474
475 // Zero out the stress tensor
476 stressTensor *= 0.0;
477 // Zero out the heatFlux
478 fDecomp_->setHeatFlux( Vector3d(0.0) );
479 }
480
481 void ForceManager::shortRangeInteractions() {
482 Molecule* mol;
483 RigidBody* rb;
484 Bond* bond;
485 Bend* bend;
486 Torsion* torsion;
487 Inversion* inversion;
488 SimInfo::MoleculeIterator mi;
489 Molecule::RigidBodyIterator rbIter;
490 Molecule::BondIterator bondIter;;
491 Molecule::BendIterator bendIter;
492 Molecule::TorsionIterator torsionIter;
493 Molecule::InversionIterator inversionIter;
494 RealType bondPotential = 0.0;
495 RealType bendPotential = 0.0;
496 RealType torsionPotential = 0.0;
497 RealType inversionPotential = 0.0;
498
499 //calculate short range interactions
500 for (mol = info_->beginMolecule(mi); mol != NULL;
501 mol = info_->nextMolecule(mi)) {
502
503 //change the positions of atoms which belong to the rigidbodies
504 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
505 rb = mol->nextRigidBody(rbIter)) {
506 rb->updateAtoms();
507 }
508
509 for (bond = mol->beginBond(bondIter); bond != NULL;
510 bond = mol->nextBond(bondIter)) {
511 bond->calcForce(doParticlePot_);
512 bondPotential += bond->getPotential();
513 }
514
515 for (bend = mol->beginBend(bendIter); bend != NULL;
516 bend = mol->nextBend(bendIter)) {
517
518 RealType angle;
519 bend->calcForce(angle, doParticlePot_);
520 RealType currBendPot = bend->getPotential();
521
522 bendPotential += bend->getPotential();
523 map<Bend*, BendDataSet>::iterator i = bendDataSets.find(bend);
524 if (i == bendDataSets.end()) {
525 BendDataSet dataSet;
526 dataSet.prev.angle = dataSet.curr.angle = angle;
527 dataSet.prev.potential = dataSet.curr.potential = currBendPot;
528 dataSet.deltaV = 0.0;
529 bendDataSets.insert(map<Bend*, BendDataSet>::value_type(bend,
530 dataSet));
531 }else {
532 i->second.prev.angle = i->second.curr.angle;
533 i->second.prev.potential = i->second.curr.potential;
534 i->second.curr.angle = angle;
535 i->second.curr.potential = currBendPot;
536 i->second.deltaV = fabs(i->second.curr.potential -
537 i->second.prev.potential);
538 }
539 }
540
541 for (torsion = mol->beginTorsion(torsionIter); torsion != NULL;
542 torsion = mol->nextTorsion(torsionIter)) {
543 RealType angle;
544 torsion->calcForce(angle, doParticlePot_);
545 RealType currTorsionPot = torsion->getPotential();
546 torsionPotential += torsion->getPotential();
547 map<Torsion*, TorsionDataSet>::iterator i = torsionDataSets.find(torsion);
548 if (i == torsionDataSets.end()) {
549 TorsionDataSet dataSet;
550 dataSet.prev.angle = dataSet.curr.angle = angle;
551 dataSet.prev.potential = dataSet.curr.potential = currTorsionPot;
552 dataSet.deltaV = 0.0;
553 torsionDataSets.insert(map<Torsion*, TorsionDataSet>::value_type(torsion, dataSet));
554 }else {
555 i->second.prev.angle = i->second.curr.angle;
556 i->second.prev.potential = i->second.curr.potential;
557 i->second.curr.angle = angle;
558 i->second.curr.potential = currTorsionPot;
559 i->second.deltaV = fabs(i->second.curr.potential -
560 i->second.prev.potential);
561 }
562 }
563
564 for (inversion = mol->beginInversion(inversionIter);
565 inversion != NULL;
566 inversion = mol->nextInversion(inversionIter)) {
567 RealType angle;
568 inversion->calcForce(angle, doParticlePot_);
569 RealType currInversionPot = inversion->getPotential();
570 inversionPotential += inversion->getPotential();
571 map<Inversion*, InversionDataSet>::iterator i = inversionDataSets.find(inversion);
572 if (i == inversionDataSets.end()) {
573 InversionDataSet dataSet;
574 dataSet.prev.angle = dataSet.curr.angle = angle;
575 dataSet.prev.potential = dataSet.curr.potential = currInversionPot;
576 dataSet.deltaV = 0.0;
577 inversionDataSets.insert(map<Inversion*, InversionDataSet>::value_type(inversion, 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 = currInversionPot;
583 i->second.deltaV = fabs(i->second.curr.potential -
584 i->second.prev.potential);
585 }
586 }
587 }
588
589 RealType shortRangePotential = bondPotential + bendPotential +
590 torsionPotential + inversionPotential;
591 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
592 curSnapshot->statData[Stats::SHORT_RANGE_POTENTIAL] = shortRangePotential;
593 curSnapshot->statData[Stats::BOND_POTENTIAL] = bondPotential;
594 curSnapshot->statData[Stats::BEND_POTENTIAL] = bendPotential;
595 curSnapshot->statData[Stats::DIHEDRAL_POTENTIAL] = torsionPotential;
596 curSnapshot->statData[Stats::INVERSION_POTENTIAL] = inversionPotential;
597 }
598
599 void ForceManager::longRangeInteractions() {
600
601
602 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
603 DataStorage* config = &(curSnapshot->atomData);
604 DataStorage* cgConfig = &(curSnapshot->cgData);
605
606 //calculate the center of mass of cutoff group
607
608 SimInfo::MoleculeIterator mi;
609 Molecule* mol;
610 Molecule::CutoffGroupIterator ci;
611 CutoffGroup* cg;
612
613 if(info_->getNCutoffGroups() > 0){
614 for (mol = info_->beginMolecule(mi); mol != NULL;
615 mol = info_->nextMolecule(mi)) {
616 for(cg = mol->beginCutoffGroup(ci); cg != NULL;
617 cg = mol->nextCutoffGroup(ci)) {
618 cg->updateCOM();
619 }
620 }
621 } else {
622 // center of mass of the group is the same as position of the atom
623 // if cutoff group does not exist
624 cgConfig->position = config->position;
625 cgConfig->velocity = config->velocity;
626 }
627
628 fDecomp_->zeroWorkArrays();
629 fDecomp_->distributeData();
630
631 int cg1, cg2, atom1, atom2, topoDist;
632 Vector3d d_grp, dag, d, gvel2, vel2;
633 RealType rgrpsq, rgrp, r2, r;
634 RealType electroMult, vdwMult;
635 RealType vij;
636 Vector3d fij, fg, f1;
637 tuple3<RealType, RealType, RealType> cuts;
638 RealType rCutSq;
639 bool in_switching_region;
640 RealType sw, dswdr, swderiv;
641 vector<int> atomListColumn, atomListRow, atomListLocal;
642 InteractionData idat;
643 SelfData sdat;
644 RealType mf;
645 RealType lrPot;
646 RealType vpair;
647 RealType dVdFQ1(0.0);
648 RealType dVdFQ2(0.0);
649 potVec longRangePotential(0.0);
650 potVec workPot(0.0);
651 vector<int>::iterator ia, jb;
652
653 int loopStart, loopEnd;
654
655 idat.vdwMult = &vdwMult;
656 idat.electroMult = &electroMult;
657 idat.pot = &workPot;
658 sdat.pot = fDecomp_->getEmbeddingPotential();
659 idat.vpair = &vpair;
660 idat.dVdFQ1 = &dVdFQ1;
661 idat.dVdFQ2 = &dVdFQ2;
662 idat.f1 = &f1;
663 idat.sw = &sw;
664 idat.shiftedPot = (cutoffMethod_ == SHIFTED_POTENTIAL) ? true : false;
665 idat.shiftedForce = (cutoffMethod_ == SHIFTED_FORCE) ? true : false;
666 idat.doParticlePot = doParticlePot_;
667 sdat.doParticlePot = doParticlePot_;
668
669 loopEnd = PAIR_LOOP;
670 if (info_->requiresPrepair() ) {
671 loopStart = PREPAIR_LOOP;
672 } else {
673 loopStart = PAIR_LOOP;
674 }
675
676 for (int iLoop = loopStart; iLoop <= loopEnd; iLoop++) {
677
678 if (iLoop == loopStart) {
679 bool update_nlist = fDecomp_->checkNeighborList();
680 if (update_nlist)
681 neighborList = fDecomp_->buildNeighborList();
682 }
683
684 for (vector<pair<int, int> >::iterator it = neighborList.begin();
685 it != neighborList.end(); ++it) {
686
687 cg1 = (*it).first;
688 cg2 = (*it).second;
689
690 cuts = fDecomp_->getGroupCutoffs(cg1, cg2);
691
692 d_grp = fDecomp_->getIntergroupVector(cg1, cg2);
693
694 curSnapshot->wrapVector(d_grp);
695 rgrpsq = d_grp.lengthSquare();
696 rCutSq = cuts.second;
697
698 if (rgrpsq < rCutSq) {
699 idat.rcut = &cuts.first;
700 if (iLoop == PAIR_LOOP) {
701 vij = 0.0;
702 fij = V3Zero;
703 }
704
705 in_switching_region = switcher_->getSwitch(rgrpsq, sw, dswdr,
706 rgrp);
707
708 atomListRow = fDecomp_->getAtomsInGroupRow(cg1);
709 atomListColumn = fDecomp_->getAtomsInGroupColumn(cg2);
710
711 if (doHeatFlux_)
712 gvel2 = fDecomp_->getGroupVelocityColumn(cg2);
713
714 for (ia = atomListRow.begin();
715 ia != atomListRow.end(); ++ia) {
716 atom1 = (*ia);
717
718 for (jb = atomListColumn.begin();
719 jb != atomListColumn.end(); ++jb) {
720 atom2 = (*jb);
721
722 if (!fDecomp_->skipAtomPair(atom1, atom2)) {
723 vpair = 0.0;
724 workPot = 0.0;
725 f1 = V3Zero;
726 dVdFQ1 = 0.0;
727 dVdFQ2 = 0.0;
728
729 fDecomp_->fillInteractionData(idat, atom1, atom2);
730
731 topoDist = fDecomp_->getTopologicalDistance(atom1, atom2);
732 vdwMult = vdwScale_[topoDist];
733 electroMult = electrostaticScale_[topoDist];
734
735 if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
736 idat.d = &d_grp;
737 idat.r2 = &rgrpsq;
738 if (doHeatFlux_)
739 vel2 = gvel2;
740 } else {
741 d = fDecomp_->getInteratomicVector(atom1, atom2);
742 curSnapshot->wrapVector( d );
743 r2 = d.lengthSquare();
744 idat.d = &d;
745 idat.r2 = &r2;
746 if (doHeatFlux_)
747 vel2 = fDecomp_->getAtomVelocityColumn(atom2);
748 }
749
750 r = sqrt( *(idat.r2) );
751 idat.rij = &r;
752
753 if (iLoop == PREPAIR_LOOP) {
754 interactionMan_->doPrePair(idat);
755 } else {
756 interactionMan_->doPair(idat);
757 fDecomp_->unpackInteractionData(idat, atom1, atom2);
758 vij += vpair;
759 fij += f1;
760 stressTensor -= outProduct( *(idat.d), f1);
761 if (doHeatFlux_)
762 fDecomp_->addToHeatFlux(*(idat.d) * dot(f1, vel2));
763 }
764 }
765 }
766 }
767
768 if (iLoop == PAIR_LOOP) {
769 if (in_switching_region) {
770 swderiv = vij * dswdr / rgrp;
771 fg = swderiv * d_grp;
772 fij += fg;
773
774 if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
775 stressTensor -= outProduct( *(idat.d), fg);
776 if (doHeatFlux_)
777 fDecomp_->addToHeatFlux(*(idat.d) * dot(fg, vel2));
778
779 }
780
781 for (ia = atomListRow.begin();
782 ia != atomListRow.end(); ++ia) {
783 atom1 = (*ia);
784 mf = fDecomp_->getMassFactorRow(atom1);
785 // fg is the force on atom ia due to cutoff group's
786 // presence in switching region
787 fg = swderiv * d_grp * mf;
788 fDecomp_->addForceToAtomRow(atom1, fg);
789 if (atomListRow.size() > 1) {
790 if (info_->usesAtomicVirial()) {
791 // find the distance between the atom
792 // and the center of the cutoff group:
793 dag = fDecomp_->getAtomToGroupVectorRow(atom1, cg1);
794 stressTensor -= outProduct(dag, fg);
795 if (doHeatFlux_)
796 fDecomp_->addToHeatFlux( dag * dot(fg, vel2));
797 }
798 }
799 }
800 for (jb = atomListColumn.begin();
801 jb != atomListColumn.end(); ++jb) {
802 atom2 = (*jb);
803 mf = fDecomp_->getMassFactorColumn(atom2);
804 // fg is the force on atom jb due to cutoff group's
805 // presence in switching region
806 fg = -swderiv * d_grp * mf;
807 fDecomp_->addForceToAtomColumn(atom2, fg);
808
809 if (atomListColumn.size() > 1) {
810 if (info_->usesAtomicVirial()) {
811 // find the distance between the atom
812 // and the center of the cutoff group:
813 dag = fDecomp_->getAtomToGroupVectorColumn(atom2, cg2);
814 stressTensor -= outProduct(dag, fg);
815 if (doHeatFlux_)
816 fDecomp_->addToHeatFlux( dag * dot(fg, vel2));
817 }
818 }
819 }
820 }
821 //if (!info_->usesAtomicVirial()) {
822 // stressTensor -= outProduct(d_grp, fij);
823 // if (doHeatFlux_)
824 // fDecomp_->addToHeatFlux( d_grp * dot(fij, vel2));
825 //}
826 }
827 }
828 }
829
830 if (iLoop == PREPAIR_LOOP) {
831 if (info_->requiresPrepair()) {
832
833 fDecomp_->collectIntermediateData();
834
835 for (int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
836 fDecomp_->fillSelfData(sdat, atom1);
837 interactionMan_->doPreForce(sdat);
838 }
839
840 fDecomp_->distributeIntermediateData();
841
842 }
843 }
844 }
845
846 fDecomp_->collectData();
847
848 if (info_->requiresSelfCorrection()) {
849
850 for (int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
851 fDecomp_->fillSelfData(sdat, atom1);
852 interactionMan_->doSelfCorrection(sdat);
853 }
854
855 }
856
857 longRangePotential = *(fDecomp_->getEmbeddingPotential()) +
858 *(fDecomp_->getPairwisePotential());
859
860 lrPot = longRangePotential.sum();
861
862 //store the stressTensor and long range potential
863 curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL] = lrPot;
864 curSnapshot->statData[Stats::VANDERWAALS_POTENTIAL] = longRangePotential[VANDERWAALS_FAMILY];
865 curSnapshot->statData[Stats::ELECTROSTATIC_POTENTIAL] = longRangePotential[ELECTROSTATIC_FAMILY];
866 }
867
868
869 void ForceManager::postCalculation() {
870 SimInfo::MoleculeIterator mi;
871 Molecule* mol;
872 Molecule::RigidBodyIterator rbIter;
873 RigidBody* rb;
874 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
875
876 // collect the atomic forces onto rigid bodies
877
878 for (mol = info_->beginMolecule(mi); mol != NULL;
879 mol = info_->nextMolecule(mi)) {
880 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
881 rb = mol->nextRigidBody(rbIter)) {
882 Mat3x3d rbTau = rb->calcForcesAndTorquesAndVirial();
883 stressTensor += rbTau;
884 }
885 }
886
887 #ifdef IS_MPI
888
889 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, stressTensor.getArrayPointer(), 9,
890 MPI::REALTYPE, MPI::SUM);
891 #endif
892 curSnapshot->setStressTensor(stressTensor);
893
894 }
895
896 } //end namespace OpenMD

Properties

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