ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/ForceManager.cpp
Revision: 1616
Committed: Tue Aug 30 15:45:35 2011 UTC (13 years, 8 months ago) by gezelter
File size: 31296 byte(s)
Log Message:
Fixing cutoff groups

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

Properties

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