ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/ForceManager.cpp
Revision: 1577
Committed: Wed Jun 8 20:26:56 2011 UTC (13 years, 10 months ago) by gezelter
File size: 26629 byte(s)
Log Message:
bug fixes

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

Properties

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