ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/ForceManager.cpp
Revision: 1613
Committed: Thu Aug 18 20:18:19 2011 UTC (13 years, 8 months ago) by gezelter
File size: 28102 byte(s)
Log Message:
Fixed a parallel bug in computing exclude lists.
Added file versioning information in MD files.
Still tracking down cutoff group bugs.

File Contents

# User Rev Content
1 gezelter 507 /*
2 gezelter 246 * 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 gezelter 1390 * 1. Redistributions of source code must retain the above copyright
10 gezelter 246 * notice, this list of conditions and the following disclaimer.
11     *
12 gezelter 1390 * 2. Redistributions in binary form must reproduce the above copyright
13 gezelter 246 * 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 gezelter 1390 *
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 gezelter 246 */
41    
42 gezelter 507 /**
43     * @file ForceManager.cpp
44     * @author tlin
45     * @date 11/09/2004
46     * @time 10:39am
47     * @version 1.0
48     */
49 gezelter 246
50 gezelter 1576
51 gezelter 246 #include "brains/ForceManager.hpp"
52     #include "primitives/Molecule.hpp"
53 gezelter 1390 #define __OPENMD_C
54 gezelter 246 #include "utils/simError.h"
55 xsun 1215 #include "primitives/Bond.hpp"
56 tim 749 #include "primitives/Bend.hpp"
57 cli2 1275 #include "primitives/Torsion.hpp"
58     #include "primitives/Inversion.hpp"
59 gezelter 1551 #include "nonbonded/NonBondedInteraction.hpp"
60 gezelter 1549 #include "parallel/ForceMatrixDecomposition.hpp"
61 gezelter 1467
62 gezelter 1583 #include <cstdio>
63     #include <iostream>
64     #include <iomanip>
65    
66 gezelter 1545 using namespace std;
67 gezelter 1390 namespace OpenMD {
68 gezelter 1469
69 gezelter 1545 ForceManager::ForceManager(SimInfo * info) : info_(info) {
70 gezelter 1576 forceField_ = info_->getForceField();
71 gezelter 1577 interactionMan_ = new InteractionManager();
72 gezelter 1579 fDecomp_ = new ForceMatrixDecomposition(info_, interactionMan_);
73 gezelter 1469 }
74 gezelter 1576
75     /**
76     * setupCutoffs
77     *
78 gezelter 1587 * Sets the values of cutoffRadius, switchingRadius, cutoffMethod,
79     * and cutoffPolicy
80 gezelter 1576 *
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 gezelter 1590 * cutoffMethod : (one of HARD, SWITCHED, SHIFTED_FORCE,
90     * or SHIFTED_POTENTIAL)
91 gezelter 1576 * 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 gezelter 1587 *
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 gezelter 1576 */
107     void ForceManager::setupCutoffs() {
108    
109     Globals* simParams_ = info_->getSimParams();
110     ForceFieldOptions& forceFieldOptions_ = forceField_->getForceFieldOptions();
111 gezelter 1613 int mdFileVersion;
112 gezelter 1576
113 gezelter 1613 if (simParams_->haveMDfileVersion())
114     mdFileVersion = simParams_->getMDfileVersion();
115     else
116     mdFileVersion = 0;
117    
118    
119 gezelter 1576 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 gezelter 1579 }
148 gezelter 1576 }
149    
150 gezelter 1583 fDecomp_->setUserCutoff(rCut_);
151 gezelter 1584 interactionMan_->setCutoffRadius(rCut_);
152 gezelter 1583
153 gezelter 1576 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 gezelter 1545
159 gezelter 1576 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     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     }
184    
185     map<string, CutoffPolicy> stringToCutoffPolicy;
186     stringToCutoffPolicy["MIX"] = MIX;
187     stringToCutoffPolicy["MAX"] = MAX;
188     stringToCutoffPolicy["TRADITIONAL"] = TRADITIONAL;
189    
190     std::string cutPolicy;
191     if (forceFieldOptions_.haveCutoffPolicy()){
192     cutPolicy = forceFieldOptions_.getCutoffPolicy();
193     }else if (simParams_->haveCutoffPolicy()) {
194     cutPolicy = simParams_->getCutoffPolicy();
195     }
196    
197     if (!cutPolicy.empty()){
198     toUpper(cutPolicy);
199     map<string, CutoffPolicy>::iterator i;
200     i = stringToCutoffPolicy.find(cutPolicy);
201    
202     if (i == stringToCutoffPolicy.end()) {
203     sprintf(painCave.errMsg,
204     "ForceManager::setupCutoffs: Could not find chosen cutoffPolicy %s\n"
205     "\tShould be one of: "
206     "MIX, MAX, or TRADITIONAL\n",
207     cutPolicy.c_str());
208     painCave.isFatal = 1;
209     painCave.severity = OPENMD_ERROR;
210     simError();
211     } else {
212     cutoffPolicy_ = i->second;
213     }
214     } else {
215     sprintf(painCave.errMsg,
216     "ForceManager::setupCutoffs: No value was set for the cutoffPolicy.\n"
217     "\tOpenMD will use TRADITIONAL.\n");
218     painCave.isFatal = 0;
219     painCave.severity = OPENMD_INFO;
220     simError();
221     cutoffPolicy_ = TRADITIONAL;
222     }
223 gezelter 1587
224 gezelter 1579 fDecomp_->setCutoffPolicy(cutoffPolicy_);
225 gezelter 1587
226     // create the switching function object:
227 gezelter 1576
228 gezelter 1577 switcher_ = new SwitchingFunction();
229 gezelter 1587
230     if (cutoffMethod_ == SWITCHED) {
231     if (simParams_->haveSwitchingRadius()) {
232     rSwitch_ = simParams_->getSwitchingRadius();
233     if (rSwitch_ > rCut_) {
234     sprintf(painCave.errMsg,
235     "ForceManager::setupCutoffs: switchingRadius (%f) is larger "
236     "than the cutoffRadius(%f)\n", rSwitch_, rCut_);
237     painCave.isFatal = 1;
238     painCave.severity = OPENMD_ERROR;
239     simError();
240     }
241     } else {
242     rSwitch_ = 0.85 * rCut_;
243 gezelter 1576 sprintf(painCave.errMsg,
244 gezelter 1587 "ForceManager::setupCutoffs: No value was set for the switchingRadius.\n"
245     "\tOpenMD will use a default value of 85 percent of the cutoffRadius.\n"
246     "\tswitchingRadius = %f. for this simulation\n", rSwitch_);
247     painCave.isFatal = 0;
248     painCave.severity = OPENMD_WARNING;
249 gezelter 1576 simError();
250     }
251 gezelter 1587 } else {
252     if (simParams_->haveSwitchingRadius()) {
253     map<string, CutoffMethod>::const_iterator it;
254     string theMeth;
255     for (it = stringToCutoffMethod.begin();
256     it != stringToCutoffMethod.end(); ++it) {
257     if (it->second == cutoffMethod_) {
258     theMeth = it->first;
259     break;
260     }
261     }
262     sprintf(painCave.errMsg,
263     "ForceManager::setupCutoffs: the cutoffMethod (%s)\n"
264     "\tis not set to SWITCHED, so switchingRadius value\n"
265     "\twill be ignored for this simulation\n", theMeth.c_str());
266     painCave.isFatal = 0;
267     painCave.severity = OPENMD_WARNING;
268     simError();
269     }
270    
271     rSwitch_ = rCut_;
272     }
273 gezelter 1576
274 gezelter 1577 // Default to cubic switching function.
275     sft_ = cubic;
276 gezelter 1576 if (simParams_->haveSwitchingFunctionType()) {
277     string funcType = simParams_->getSwitchingFunctionType();
278     toUpper(funcType);
279     if (funcType == "CUBIC") {
280     sft_ = cubic;
281     } else {
282     if (funcType == "FIFTH_ORDER_POLYNOMIAL") {
283     sft_ = fifth_order_poly;
284     } else {
285     // throw error
286     sprintf( painCave.errMsg,
287     "ForceManager::setupSwitching : Unknown switchingFunctionType. (Input file specified %s .)\n"
288     "\tswitchingFunctionType must be one of: "
289     "\"cubic\" or \"fifth_order_polynomial\".",
290     funcType.c_str() );
291     painCave.isFatal = 1;
292     painCave.severity = OPENMD_ERROR;
293     simError();
294     }
295     }
296     }
297     switcher_->setSwitchType(sft_);
298     switcher_->setSwitch(rSwitch_, rCut_);
299 gezelter 1584 interactionMan_->setSwitchingRadius(rSwitch_);
300 gezelter 1576 }
301    
302     void ForceManager::initialize() {
303    
304 gezelter 1569 if (!info_->isTopologyDone()) {
305 gezelter 1590
306 gezelter 507 info_->update();
307 gezelter 1546 interactionMan_->setSimInfo(info_);
308     interactionMan_->initialize();
309 gezelter 1576
310     // We want to delay the cutoffs until after the interaction
311     // manager has set up the atom-atom interactions so that we can
312     // query them for suggested cutoff values
313     setupCutoffs();
314    
315     info_->prepareTopology();
316 gezelter 246 }
317 gezelter 1576
318     ForceFieldOptions& fopts = forceField_->getForceFieldOptions();
319 gezelter 1126
320 gezelter 1590 // Force fields can set options on how to scale van der Waals and
321     // electrostatic interactions for atoms connected via bonds, bends
322     // and torsions in this case the topological distance between
323     // atoms is:
324 gezelter 1576 // 0 = topologically unconnected
325     // 1 = bonded together
326     // 2 = connected via a bend
327     // 3 = connected via a torsion
328    
329     vdwScale_.reserve(4);
330     fill(vdwScale_.begin(), vdwScale_.end(), 0.0);
331    
332     electrostaticScale_.reserve(4);
333     fill(electrostaticScale_.begin(), electrostaticScale_.end(), 0.0);
334    
335     vdwScale_[0] = 1.0;
336     vdwScale_[1] = fopts.getvdw12scale();
337     vdwScale_[2] = fopts.getvdw13scale();
338     vdwScale_[3] = fopts.getvdw14scale();
339    
340     electrostaticScale_[0] = 1.0;
341     electrostaticScale_[1] = fopts.getelectrostatic12scale();
342     electrostaticScale_[2] = fopts.getelectrostatic13scale();
343     electrostaticScale_[3] = fopts.getelectrostatic14scale();
344    
345     fDecomp_->distributeInitialData();
346    
347     initialized_ = true;
348    
349     }
350    
351     void ForceManager::calcForces() {
352    
353     if (!initialized_) initialize();
354    
355 gezelter 1544 preCalculation();
356 gezelter 1546 shortRangeInteractions();
357     longRangeInteractions();
358 gezelter 1576 postCalculation();
359 gezelter 507 }
360 gezelter 1126
361 gezelter 507 void ForceManager::preCalculation() {
362 gezelter 246 SimInfo::MoleculeIterator mi;
363     Molecule* mol;
364     Molecule::AtomIterator ai;
365     Atom* atom;
366     Molecule::RigidBodyIterator rbIter;
367     RigidBody* rb;
368 gezelter 1540 Molecule::CutoffGroupIterator ci;
369     CutoffGroup* cg;
370 gezelter 246
371     // forces are zeroed here, before any are accumulated.
372 chuckv 1245
373 gezelter 1126 for (mol = info_->beginMolecule(mi); mol != NULL;
374     mol = info_->nextMolecule(mi)) {
375 gezelter 1590 for(atom = mol->beginAtom(ai); atom != NULL;
376     atom = mol->nextAtom(ai)) {
377 gezelter 507 atom->zeroForcesAndTorques();
378     }
379 gezelter 1590
380 gezelter 507 //change the positions of atoms which belong to the rigidbodies
381 gezelter 1126 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
382     rb = mol->nextRigidBody(rbIter)) {
383 gezelter 507 rb->zeroForcesAndTorques();
384     }
385 gezelter 1590
386 gezelter 1540 if(info_->getNGlobalCutoffGroups() != info_->getNGlobalAtoms()){
387     for(cg = mol->beginCutoffGroup(ci); cg != NULL;
388     cg = mol->nextCutoffGroup(ci)) {
389     //calculate the center of mass of cutoff group
390     cg->updateCOM();
391     }
392     }
393 gezelter 246 }
394 gezelter 1590
395 gezelter 1126 // Zero out the stress tensor
396 gezelter 1591 tau *= 0.0;
397 gezelter 1126
398 gezelter 507 }
399 gezelter 1126
400 gezelter 1546 void ForceManager::shortRangeInteractions() {
401 gezelter 246 Molecule* mol;
402     RigidBody* rb;
403     Bond* bond;
404     Bend* bend;
405     Torsion* torsion;
406 cli2 1275 Inversion* inversion;
407 gezelter 246 SimInfo::MoleculeIterator mi;
408     Molecule::RigidBodyIterator rbIter;
409     Molecule::BondIterator bondIter;;
410     Molecule::BendIterator bendIter;
411     Molecule::TorsionIterator torsionIter;
412 cli2 1275 Molecule::InversionIterator inversionIter;
413 tim 963 RealType bondPotential = 0.0;
414     RealType bendPotential = 0.0;
415     RealType torsionPotential = 0.0;
416 cli2 1275 RealType inversionPotential = 0.0;
417 gezelter 246
418     //calculate short range interactions
419 gezelter 1126 for (mol = info_->beginMolecule(mi); mol != NULL;
420     mol = info_->nextMolecule(mi)) {
421 gezelter 246
422 gezelter 507 //change the positions of atoms which belong to the rigidbodies
423 gezelter 1126 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
424     rb = mol->nextRigidBody(rbIter)) {
425     rb->updateAtoms();
426 gezelter 507 }
427 gezelter 246
428 gezelter 1126 for (bond = mol->beginBond(bondIter); bond != NULL;
429     bond = mol->nextBond(bondIter)) {
430 tim 749 bond->calcForce();
431     bondPotential += bond->getPotential();
432 gezelter 507 }
433 gezelter 246
434 gezelter 1126 for (bend = mol->beginBend(bendIter); bend != NULL;
435     bend = mol->nextBend(bendIter)) {
436    
437     RealType angle;
438     bend->calcForce(angle);
439     RealType currBendPot = bend->getPotential();
440 gezelter 1448
441 gezelter 1126 bendPotential += bend->getPotential();
442 gezelter 1545 map<Bend*, BendDataSet>::iterator i = bendDataSets.find(bend);
443 gezelter 1126 if (i == bendDataSets.end()) {
444     BendDataSet dataSet;
445     dataSet.prev.angle = dataSet.curr.angle = angle;
446     dataSet.prev.potential = dataSet.curr.potential = currBendPot;
447     dataSet.deltaV = 0.0;
448 gezelter 1590 bendDataSets.insert(map<Bend*, BendDataSet>::value_type(bend,
449     dataSet));
450 gezelter 1126 }else {
451     i->second.prev.angle = i->second.curr.angle;
452     i->second.prev.potential = i->second.curr.potential;
453     i->second.curr.angle = angle;
454     i->second.curr.potential = currBendPot;
455     i->second.deltaV = fabs(i->second.curr.potential -
456     i->second.prev.potential);
457     }
458 gezelter 507 }
459 gezelter 1126
460     for (torsion = mol->beginTorsion(torsionIter); torsion != NULL;
461     torsion = mol->nextTorsion(torsionIter)) {
462 tim 963 RealType angle;
463 gezelter 1126 torsion->calcForce(angle);
464 tim 963 RealType currTorsionPot = torsion->getPotential();
465 gezelter 1126 torsionPotential += torsion->getPotential();
466 gezelter 1545 map<Torsion*, TorsionDataSet>::iterator i = torsionDataSets.find(torsion);
467 gezelter 1126 if (i == torsionDataSets.end()) {
468     TorsionDataSet dataSet;
469     dataSet.prev.angle = dataSet.curr.angle = angle;
470     dataSet.prev.potential = dataSet.curr.potential = currTorsionPot;
471     dataSet.deltaV = 0.0;
472 gezelter 1545 torsionDataSets.insert(map<Torsion*, TorsionDataSet>::value_type(torsion, dataSet));
473 gezelter 1126 }else {
474     i->second.prev.angle = i->second.curr.angle;
475     i->second.prev.potential = i->second.curr.potential;
476     i->second.curr.angle = angle;
477     i->second.curr.potential = currTorsionPot;
478     i->second.deltaV = fabs(i->second.curr.potential -
479     i->second.prev.potential);
480     }
481     }
482 gezelter 1545
483 cli2 1275 for (inversion = mol->beginInversion(inversionIter);
484     inversion != NULL;
485     inversion = mol->nextInversion(inversionIter)) {
486     RealType angle;
487     inversion->calcForce(angle);
488     RealType currInversionPot = inversion->getPotential();
489     inversionPotential += inversion->getPotential();
490 gezelter 1545 map<Inversion*, InversionDataSet>::iterator i = inversionDataSets.find(inversion);
491 cli2 1275 if (i == inversionDataSets.end()) {
492     InversionDataSet dataSet;
493     dataSet.prev.angle = dataSet.curr.angle = angle;
494     dataSet.prev.potential = dataSet.curr.potential = currInversionPot;
495     dataSet.deltaV = 0.0;
496 gezelter 1545 inversionDataSets.insert(map<Inversion*, InversionDataSet>::value_type(inversion, dataSet));
497 cli2 1275 }else {
498     i->second.prev.angle = i->second.curr.angle;
499     i->second.prev.potential = i->second.curr.potential;
500     i->second.curr.angle = angle;
501     i->second.curr.potential = currInversionPot;
502     i->second.deltaV = fabs(i->second.curr.potential -
503     i->second.prev.potential);
504     }
505     }
506 gezelter 246 }
507    
508 gezelter 1126 RealType shortRangePotential = bondPotential + bendPotential +
509 cli2 1275 torsionPotential + inversionPotential;
510 gezelter 246 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
511     curSnapshot->statData[Stats::SHORT_RANGE_POTENTIAL] = shortRangePotential;
512 tim 665 curSnapshot->statData[Stats::BOND_POTENTIAL] = bondPotential;
513     curSnapshot->statData[Stats::BEND_POTENTIAL] = bendPotential;
514     curSnapshot->statData[Stats::DIHEDRAL_POTENTIAL] = torsionPotential;
515 gezelter 1545 curSnapshot->statData[Stats::INVERSION_POTENTIAL] = inversionPotential;
516 gezelter 507 }
517 gezelter 1126
518 gezelter 1546 void ForceManager::longRangeInteractions() {
519 gezelter 1581
520 gezelter 1545 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
521     DataStorage* config = &(curSnapshot->atomData);
522     DataStorage* cgConfig = &(curSnapshot->cgData);
523    
524 gezelter 1581 //calculate the center of mass of cutoff group
525    
526     SimInfo::MoleculeIterator mi;
527     Molecule* mol;
528     Molecule::CutoffGroupIterator ci;
529     CutoffGroup* cg;
530    
531     if(info_->getNCutoffGroups() > 0){
532     for (mol = info_->beginMolecule(mi); mol != NULL;
533     mol = info_->nextMolecule(mi)) {
534     for(cg = mol->beginCutoffGroup(ci); cg != NULL;
535     cg = mol->nextCutoffGroup(ci)) {
536     cg->updateCOM();
537     }
538     }
539     } else {
540     // center of mass of the group is the same as position of the atom
541     // if cutoff group does not exist
542     cgConfig->position = config->position;
543     }
544    
545 gezelter 1575 fDecomp_->zeroWorkArrays();
546 gezelter 1549 fDecomp_->distributeData();
547 gezelter 1579
548     int cg1, cg2, atom1, atom2, topoDist;
549     Vector3d d_grp, dag, d;
550     RealType rgrpsq, rgrp, r2, r;
551     RealType electroMult, vdwMult;
552 gezelter 1549 RealType vij;
553 gezelter 1581 Vector3d fij, fg, f1;
554 gezelter 1576 tuple3<RealType, RealType, RealType> cuts;
555 gezelter 1545 RealType rCutSq;
556     bool in_switching_region;
557     RealType sw, dswdr, swderiv;
558 gezelter 1549 vector<int> atomListColumn, atomListRow, atomListLocal;
559 gezelter 1545 InteractionData idat;
560 gezelter 1546 SelfData sdat;
561     RealType mf;
562 gezelter 1575 RealType lrPot;
563 gezelter 1579 RealType vpair;
564 gezelter 1583 potVec longRangePotential(0.0);
565     potVec workPot(0.0);
566 gezelter 1544
567 gezelter 1545 int loopStart, loopEnd;
568 gezelter 1544
569 gezelter 1581 idat.vdwMult = &vdwMult;
570     idat.electroMult = &electroMult;
571 gezelter 1583 idat.pot = &workPot;
572     sdat.pot = fDecomp_->getEmbeddingPotential();
573 gezelter 1581 idat.vpair = &vpair;
574     idat.f1 = &f1;
575     idat.sw = &sw;
576 gezelter 1583 idat.shiftedPot = (cutoffMethod_ == SHIFTED_POTENTIAL) ? true : false;
577     idat.shiftedForce = (cutoffMethod_ == SHIFTED_FORCE) ? true : false;
578    
579 gezelter 1545 loopEnd = PAIR_LOOP;
580 gezelter 1546 if (info_->requiresPrepair() ) {
581 gezelter 1545 loopStart = PREPAIR_LOOP;
582     } else {
583     loopStart = PAIR_LOOP;
584     }
585 gezelter 1583
586 gezelter 1579 for (int iLoop = loopStart; iLoop <= loopEnd; iLoop++) {
587    
588 gezelter 1545 if (iLoop == loopStart) {
589 gezelter 1549 bool update_nlist = fDecomp_->checkNeighborList();
590 gezelter 1545 if (update_nlist)
591 gezelter 1549 neighborList = fDecomp_->buildNeighborList();
592 gezelter 1612 }
593    
594 gezelter 1545 for (vector<pair<int, int> >::iterator it = neighborList.begin();
595     it != neighborList.end(); ++it) {
596 gezelter 1579
597 gezelter 1545 cg1 = (*it).first;
598     cg2 = (*it).second;
599 gezelter 1576
600     cuts = fDecomp_->getGroupCutoffs(cg1, cg2);
601 gezelter 1545
602 gezelter 1549 d_grp = fDecomp_->getIntergroupVector(cg1, cg2);
603 gezelter 1613
604 gezelter 1545 curSnapshot->wrapVector(d_grp);
605     rgrpsq = d_grp.lengthSquare();
606 gezelter 1576 rCutSq = cuts.second;
607    
608 gezelter 1545 if (rgrpsq < rCutSq) {
609 gezelter 1579 idat.rcut = &cuts.first;
610 gezelter 1545 if (iLoop == PAIR_LOOP) {
611 gezelter 1587 vij = 0.0;
612 gezelter 1545 fij = V3Zero;
613     }
614    
615 gezelter 1579 in_switching_region = switcher_->getSwitch(rgrpsq, sw, dswdr,
616 gezelter 1576 rgrp);
617    
618 gezelter 1549 atomListRow = fDecomp_->getAtomsInGroupRow(cg1);
619     atomListColumn = fDecomp_->getAtomsInGroupColumn(cg2);
620 gezelter 1613
621 gezelter 1545
622 gezelter 1549 for (vector<int>::iterator ia = atomListRow.begin();
623     ia != atomListRow.end(); ++ia) {
624 gezelter 1545 atom1 = (*ia);
625    
626 gezelter 1549 for (vector<int>::iterator jb = atomListColumn.begin();
627     jb != atomListColumn.end(); ++jb) {
628 gezelter 1545 atom2 = (*jb);
629 gezelter 1593
630 gezelter 1549 if (!fDecomp_->skipAtomPair(atom1, atom2)) {
631 gezelter 1579 vpair = 0.0;
632 gezelter 1583 workPot = 0.0;
633 gezelter 1581 f1 = V3Zero;
634 gezelter 1575
635 gezelter 1581 fDecomp_->fillInteractionData(idat, atom1, atom2);
636 gezelter 1579
637     topoDist = fDecomp_->getTopologicalDistance(atom1, atom2);
638     vdwMult = vdwScale_[topoDist];
639     electroMult = electrostaticScale_[topoDist];
640 gezelter 1546
641 gezelter 1549 if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
642 gezelter 1579 idat.d = &d_grp;
643     idat.r2 = &rgrpsq;
644 gezelter 1545 } else {
645 gezelter 1579 d = fDecomp_->getInteratomicVector(atom1, atom2);
646     curSnapshot->wrapVector( d );
647     r2 = d.lengthSquare();
648     idat.d = &d;
649     idat.r2 = &r2;
650 gezelter 1545 }
651 gezelter 1601
652 gezelter 1581 r = sqrt( *(idat.r2) );
653 gezelter 1579 idat.rij = &r;
654 gezelter 1546
655 gezelter 1545 if (iLoop == PREPAIR_LOOP) {
656     interactionMan_->doPrePair(idat);
657     } else {
658     interactionMan_->doPair(idat);
659 gezelter 1575 fDecomp_->unpackInteractionData(idat, atom1, atom2);
660 gezelter 1593
661 gezelter 1581 vij += vpair;
662     fij += f1;
663     tau -= outProduct( *(idat.d), f1);
664 gezelter 1545 }
665     }
666     }
667     }
668    
669     if (iLoop == PAIR_LOOP) {
670     if (in_switching_region) {
671     swderiv = vij * dswdr / rgrp;
672     fg = swderiv * d_grp;
673     fij += fg;
674    
675 gezelter 1549 if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
676 gezelter 1554 tau -= outProduct( *(idat.d), fg);
677 gezelter 1545 }
678    
679 gezelter 1549 for (vector<int>::iterator ia = atomListRow.begin();
680     ia != atomListRow.end(); ++ia) {
681 gezelter 1545 atom1 = (*ia);
682 gezelter 1569 mf = fDecomp_->getMassFactorRow(atom1);
683 gezelter 1545 // fg is the force on atom ia due to cutoff group's
684     // presence in switching region
685     fg = swderiv * d_grp * mf;
686 gezelter 1549 fDecomp_->addForceToAtomRow(atom1, fg);
687 gezelter 1545
688 gezelter 1549 if (atomListRow.size() > 1) {
689 gezelter 1546 if (info_->usesAtomicVirial()) {
690 gezelter 1545 // find the distance between the atom
691     // and the center of the cutoff group:
692 gezelter 1549 dag = fDecomp_->getAtomToGroupVectorRow(atom1, cg1);
693 gezelter 1545 tau -= outProduct(dag, fg);
694     }
695     }
696     }
697 gezelter 1549 for (vector<int>::iterator jb = atomListColumn.begin();
698     jb != atomListColumn.end(); ++jb) {
699 gezelter 1545 atom2 = (*jb);
700 gezelter 1569 mf = fDecomp_->getMassFactorColumn(atom2);
701 gezelter 1545 // fg is the force on atom jb due to cutoff group's
702     // presence in switching region
703     fg = -swderiv * d_grp * mf;
704 gezelter 1549 fDecomp_->addForceToAtomColumn(atom2, fg);
705 gezelter 1545
706 gezelter 1549 if (atomListColumn.size() > 1) {
707 gezelter 1546 if (info_->usesAtomicVirial()) {
708 gezelter 1545 // find the distance between the atom
709     // and the center of the cutoff group:
710 gezelter 1549 dag = fDecomp_->getAtomToGroupVectorColumn(atom2, cg2);
711 gezelter 1545 tau -= outProduct(dag, fg);
712     }
713     }
714     }
715     }
716 gezelter 1613 //if (!info_->usesAtomicVirial()) {
717 gezelter 1545 // tau -= outProduct(d_grp, fij);
718     //}
719     }
720     }
721     }
722    
723     if (iLoop == PREPAIR_LOOP) {
724 gezelter 1590 if (info_->requiresPrepair()) {
725    
726 gezelter 1549 fDecomp_->collectIntermediateData();
727 gezelter 1570
728     for (int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
729 gezelter 1581 fDecomp_->fillSelfData(sdat, atom1);
730 gezelter 1545 interactionMan_->doPreForce(sdat);
731     }
732 gezelter 1590
733     fDecomp_->distributeIntermediateData();
734    
735 gezelter 1545 }
736     }
737 gezelter 1544 }
738 gezelter 1545
739 gezelter 1549 fDecomp_->collectData();
740 gezelter 1570
741     if (info_->requiresSelfCorrection()) {
742 gezelter 1545
743 gezelter 1570 for (int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
744 gezelter 1581 fDecomp_->fillSelfData(sdat, atom1);
745 gezelter 1570 interactionMan_->doSelfCorrection(sdat);
746     }
747    
748     }
749    
750 gezelter 1583 longRangePotential = *(fDecomp_->getEmbeddingPotential()) +
751     *(fDecomp_->getPairwisePotential());
752    
753 gezelter 1575 lrPot = longRangePotential.sum();
754    
755 gezelter 246 //store the tau and long range potential
756 chuckv 664 curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL] = lrPot;
757 gezelter 1550 curSnapshot->statData[Stats::VANDERWAALS_POTENTIAL] = longRangePotential[VANDERWAALS_FAMILY];
758     curSnapshot->statData[Stats::ELECTROSTATIC_POTENTIAL] = longRangePotential[ELECTROSTATIC_FAMILY];
759 gezelter 507 }
760 gezelter 246
761 gezelter 1126
762 gezelter 1464 void ForceManager::postCalculation() {
763 gezelter 246 SimInfo::MoleculeIterator mi;
764     Molecule* mol;
765     Molecule::RigidBodyIterator rbIter;
766     RigidBody* rb;
767 gezelter 1126 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
768 gezelter 246
769     // collect the atomic forces onto rigid bodies
770 gezelter 1126
771     for (mol = info_->beginMolecule(mi); mol != NULL;
772     mol = info_->nextMolecule(mi)) {
773     for (rb = mol->beginRigidBody(rbIter); rb != NULL;
774     rb = mol->nextRigidBody(rbIter)) {
775 gezelter 1464 Mat3x3d rbTau = rb->calcForcesAndTorquesAndVirial();
776     tau += rbTau;
777 gezelter 507 }
778 gezelter 1126 }
779 gezelter 1464
780 gezelter 1126 #ifdef IS_MPI
781 gezelter 1464 Mat3x3d tmpTau(tau);
782     MPI_Allreduce(tmpTau.getArrayPointer(), tau.getArrayPointer(),
783     9, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
784 gezelter 1126 #endif
785 gezelter 1464 curSnapshot->statData.setTau(tau);
786 gezelter 507 }
787 gezelter 246
788 gezelter 1390 } //end namespace OpenMD

Properties

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