ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/devel_omp/src/brains/ForceManager.cpp
Revision: 1595
Committed: Tue Jul 19 18:50:04 2011 UTC (13 years, 9 months ago) by chuckv
File size: 35589 byte(s)
Log Message:
Adding initial OpenMP support using new neighbor lists.


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    
112     if (simParams_->haveCutoffRadius()) {
113     rCut_ = simParams_->getCutoffRadius();
114     } else {
115     if (info_->usesElectrostaticAtoms()) {
116     sprintf(painCave.errMsg,
117     "ForceManager::setupCutoffs: No value was set for the cutoffRadius.\n"
118     "\tOpenMD will use a default value of 12.0 angstroms"
119     "\tfor the cutoffRadius.\n");
120     painCave.isFatal = 0;
121     painCave.severity = OPENMD_INFO;
122     simError();
123     rCut_ = 12.0;
124     } else {
125     RealType thisCut;
126     set<AtomType*>::iterator i;
127     set<AtomType*> atomTypes;
128     atomTypes = info_->getSimulatedAtomTypes();
129     for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
130     thisCut = interactionMan_->getSuggestedCutoffRadius((*i));
131     rCut_ = max(thisCut, rCut_);
132     }
133     sprintf(painCave.errMsg,
134     "ForceManager::setupCutoffs: No value was set for the cutoffRadius.\n"
135     "\tOpenMD will use %lf angstroms.\n",
136     rCut_);
137     painCave.isFatal = 0;
138     painCave.severity = OPENMD_INFO;
139     simError();
140 gezelter 1579 }
141 gezelter 1576 }
142    
143 gezelter 1583 fDecomp_->setUserCutoff(rCut_);
144 gezelter 1584 interactionMan_->setCutoffRadius(rCut_);
145 gezelter 1583
146 gezelter 1576 map<string, CutoffMethod> stringToCutoffMethod;
147     stringToCutoffMethod["HARD"] = HARD;
148     stringToCutoffMethod["SWITCHED"] = SWITCHED;
149     stringToCutoffMethod["SHIFTED_POTENTIAL"] = SHIFTED_POTENTIAL;
150     stringToCutoffMethod["SHIFTED_FORCE"] = SHIFTED_FORCE;
151 gezelter 1545
152 gezelter 1576 if (simParams_->haveCutoffMethod()) {
153     string cutMeth = toUpperCopy(simParams_->getCutoffMethod());
154     map<string, CutoffMethod>::iterator i;
155     i = stringToCutoffMethod.find(cutMeth);
156     if (i == stringToCutoffMethod.end()) {
157     sprintf(painCave.errMsg,
158     "ForceManager::setupCutoffs: Could not find chosen cutoffMethod %s\n"
159     "\tShould be one of: "
160     "HARD, SWITCHED, SHIFTED_POTENTIAL, or SHIFTED_FORCE\n",
161     cutMeth.c_str());
162     painCave.isFatal = 1;
163     painCave.severity = OPENMD_ERROR;
164     simError();
165     } else {
166     cutoffMethod_ = i->second;
167     }
168     } else {
169     sprintf(painCave.errMsg,
170     "ForceManager::setupCutoffs: No value was set for the cutoffMethod.\n"
171     "\tOpenMD will use SHIFTED_FORCE.\n");
172     painCave.isFatal = 0;
173     painCave.severity = OPENMD_INFO;
174     simError();
175     cutoffMethod_ = SHIFTED_FORCE;
176     }
177    
178     map<string, CutoffPolicy> stringToCutoffPolicy;
179     stringToCutoffPolicy["MIX"] = MIX;
180     stringToCutoffPolicy["MAX"] = MAX;
181     stringToCutoffPolicy["TRADITIONAL"] = TRADITIONAL;
182    
183     std::string cutPolicy;
184     if (forceFieldOptions_.haveCutoffPolicy()){
185     cutPolicy = forceFieldOptions_.getCutoffPolicy();
186     }else if (simParams_->haveCutoffPolicy()) {
187     cutPolicy = simParams_->getCutoffPolicy();
188     }
189    
190     if (!cutPolicy.empty()){
191     toUpper(cutPolicy);
192     map<string, CutoffPolicy>::iterator i;
193     i = stringToCutoffPolicy.find(cutPolicy);
194    
195     if (i == stringToCutoffPolicy.end()) {
196     sprintf(painCave.errMsg,
197     "ForceManager::setupCutoffs: Could not find chosen cutoffPolicy %s\n"
198     "\tShould be one of: "
199     "MIX, MAX, or TRADITIONAL\n",
200     cutPolicy.c_str());
201     painCave.isFatal = 1;
202     painCave.severity = OPENMD_ERROR;
203     simError();
204     } else {
205     cutoffPolicy_ = i->second;
206     }
207     } else {
208     sprintf(painCave.errMsg,
209     "ForceManager::setupCutoffs: No value was set for the cutoffPolicy.\n"
210     "\tOpenMD will use TRADITIONAL.\n");
211     painCave.isFatal = 0;
212     painCave.severity = OPENMD_INFO;
213     simError();
214     cutoffPolicy_ = TRADITIONAL;
215     }
216 gezelter 1587
217 gezelter 1579 fDecomp_->setCutoffPolicy(cutoffPolicy_);
218 gezelter 1587
219     // create the switching function object:
220 gezelter 1576
221 gezelter 1577 switcher_ = new SwitchingFunction();
222 gezelter 1587
223     if (cutoffMethod_ == SWITCHED) {
224     if (simParams_->haveSwitchingRadius()) {
225     rSwitch_ = simParams_->getSwitchingRadius();
226     if (rSwitch_ > rCut_) {
227     sprintf(painCave.errMsg,
228     "ForceManager::setupCutoffs: switchingRadius (%f) is larger "
229     "than the cutoffRadius(%f)\n", rSwitch_, rCut_);
230     painCave.isFatal = 1;
231     painCave.severity = OPENMD_ERROR;
232     simError();
233     }
234     } else {
235     rSwitch_ = 0.85 * rCut_;
236 gezelter 1576 sprintf(painCave.errMsg,
237 gezelter 1587 "ForceManager::setupCutoffs: No value was set for the switchingRadius.\n"
238     "\tOpenMD will use a default value of 85 percent of the cutoffRadius.\n"
239     "\tswitchingRadius = %f. for this simulation\n", rSwitch_);
240     painCave.isFatal = 0;
241     painCave.severity = OPENMD_WARNING;
242 gezelter 1576 simError();
243     }
244 gezelter 1587 } else {
245     if (simParams_->haveSwitchingRadius()) {
246     map<string, CutoffMethod>::const_iterator it;
247     string theMeth;
248     for (it = stringToCutoffMethod.begin();
249     it != stringToCutoffMethod.end(); ++it) {
250     if (it->second == cutoffMethod_) {
251     theMeth = it->first;
252     break;
253     }
254     }
255     sprintf(painCave.errMsg,
256     "ForceManager::setupCutoffs: the cutoffMethod (%s)\n"
257     "\tis not set to SWITCHED, so switchingRadius value\n"
258     "\twill be ignored for this simulation\n", theMeth.c_str());
259     painCave.isFatal = 0;
260     painCave.severity = OPENMD_WARNING;
261     simError();
262     }
263    
264     rSwitch_ = rCut_;
265     }
266 gezelter 1576
267 gezelter 1577 // Default to cubic switching function.
268     sft_ = cubic;
269 gezelter 1576 if (simParams_->haveSwitchingFunctionType()) {
270     string funcType = simParams_->getSwitchingFunctionType();
271     toUpper(funcType);
272     if (funcType == "CUBIC") {
273     sft_ = cubic;
274     } else {
275     if (funcType == "FIFTH_ORDER_POLYNOMIAL") {
276     sft_ = fifth_order_poly;
277     } else {
278     // throw error
279     sprintf( painCave.errMsg,
280     "ForceManager::setupSwitching : Unknown switchingFunctionType. (Input file specified %s .)\n"
281     "\tswitchingFunctionType must be one of: "
282     "\"cubic\" or \"fifth_order_polynomial\".",
283     funcType.c_str() );
284     painCave.isFatal = 1;
285     painCave.severity = OPENMD_ERROR;
286     simError();
287     }
288     }
289     }
290     switcher_->setSwitchType(sft_);
291     switcher_->setSwitch(rSwitch_, rCut_);
292 gezelter 1584 interactionMan_->setSwitchingRadius(rSwitch_);
293 gezelter 1576 }
294    
295     void ForceManager::initialize() {
296    
297 gezelter 1569 if (!info_->isTopologyDone()) {
298 gezelter 1590
299 gezelter 507 info_->update();
300 gezelter 1546 interactionMan_->setSimInfo(info_);
301     interactionMan_->initialize();
302 gezelter 1576
303     // We want to delay the cutoffs until after the interaction
304     // manager has set up the atom-atom interactions so that we can
305     // query them for suggested cutoff values
306     setupCutoffs();
307    
308     info_->prepareTopology();
309 gezelter 246 }
310 gezelter 1576
311     ForceFieldOptions& fopts = forceField_->getForceFieldOptions();
312 gezelter 1126
313 gezelter 1590 // Force fields can set options on how to scale van der Waals and
314     // electrostatic interactions for atoms connected via bonds, bends
315     // and torsions in this case the topological distance between
316     // atoms is:
317 gezelter 1576 // 0 = topologically unconnected
318     // 1 = bonded together
319     // 2 = connected via a bend
320     // 3 = connected via a torsion
321    
322     vdwScale_.reserve(4);
323     fill(vdwScale_.begin(), vdwScale_.end(), 0.0);
324    
325     electrostaticScale_.reserve(4);
326     fill(electrostaticScale_.begin(), electrostaticScale_.end(), 0.0);
327    
328     vdwScale_[0] = 1.0;
329     vdwScale_[1] = fopts.getvdw12scale();
330     vdwScale_[2] = fopts.getvdw13scale();
331     vdwScale_[3] = fopts.getvdw14scale();
332    
333     electrostaticScale_[0] = 1.0;
334     electrostaticScale_[1] = fopts.getelectrostatic12scale();
335     electrostaticScale_[2] = fopts.getelectrostatic13scale();
336     electrostaticScale_[3] = fopts.getelectrostatic14scale();
337    
338     fDecomp_->distributeInitialData();
339    
340     initialized_ = true;
341    
342     }
343    
344     void ForceManager::calcForces() {
345    
346     if (!initialized_) initialize();
347    
348 gezelter 1544 preCalculation();
349 gezelter 1546 shortRangeInteractions();
350 chuckv 1595 // longRangeInteractions();
351     longRangeInteractionsRapaport();
352 gezelter 1576 postCalculation();
353 gezelter 507 }
354 gezelter 1126
355 gezelter 507 void ForceManager::preCalculation() {
356 gezelter 246 SimInfo::MoleculeIterator mi;
357     Molecule* mol;
358     Molecule::AtomIterator ai;
359     Atom* atom;
360     Molecule::RigidBodyIterator rbIter;
361     RigidBody* rb;
362 gezelter 1540 Molecule::CutoffGroupIterator ci;
363     CutoffGroup* cg;
364 gezelter 246
365     // forces are zeroed here, before any are accumulated.
366 chuckv 1245
367 gezelter 1126 for (mol = info_->beginMolecule(mi); mol != NULL;
368     mol = info_->nextMolecule(mi)) {
369 gezelter 1590 for(atom = mol->beginAtom(ai); atom != NULL;
370     atom = mol->nextAtom(ai)) {
371 gezelter 507 atom->zeroForcesAndTorques();
372     }
373 gezelter 1590
374 gezelter 507 //change the positions of atoms which belong to the rigidbodies
375 gezelter 1126 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
376     rb = mol->nextRigidBody(rbIter)) {
377 gezelter 507 rb->zeroForcesAndTorques();
378     }
379 gezelter 1590
380 gezelter 1540 if(info_->getNGlobalCutoffGroups() != info_->getNGlobalAtoms()){
381     for(cg = mol->beginCutoffGroup(ci); cg != NULL;
382     cg = mol->nextCutoffGroup(ci)) {
383     //calculate the center of mass of cutoff group
384     cg->updateCOM();
385     }
386     }
387 gezelter 246 }
388 gezelter 1590
389 gezelter 1126 // Zero out the stress tensor
390 gezelter 1591 tau *= 0.0;
391 gezelter 1126
392 gezelter 507 }
393 gezelter 1126
394 gezelter 1546 void ForceManager::shortRangeInteractions() {
395 gezelter 246 Molecule* mol;
396     RigidBody* rb;
397     Bond* bond;
398     Bend* bend;
399     Torsion* torsion;
400 cli2 1275 Inversion* inversion;
401 gezelter 246 SimInfo::MoleculeIterator mi;
402     Molecule::RigidBodyIterator rbIter;
403     Molecule::BondIterator bondIter;;
404     Molecule::BendIterator bendIter;
405     Molecule::TorsionIterator torsionIter;
406 cli2 1275 Molecule::InversionIterator inversionIter;
407 tim 963 RealType bondPotential = 0.0;
408     RealType bendPotential = 0.0;
409     RealType torsionPotential = 0.0;
410 cli2 1275 RealType inversionPotential = 0.0;
411 gezelter 246
412     //calculate short range interactions
413 gezelter 1126 for (mol = info_->beginMolecule(mi); mol != NULL;
414     mol = info_->nextMolecule(mi)) {
415 gezelter 246
416 gezelter 507 //change the positions of atoms which belong to the rigidbodies
417 gezelter 1126 for (rb = mol->beginRigidBody(rbIter); rb != NULL;
418     rb = mol->nextRigidBody(rbIter)) {
419     rb->updateAtoms();
420 gezelter 507 }
421 gezelter 246
422 gezelter 1126 for (bond = mol->beginBond(bondIter); bond != NULL;
423     bond = mol->nextBond(bondIter)) {
424 tim 749 bond->calcForce();
425     bondPotential += bond->getPotential();
426 gezelter 507 }
427 gezelter 246
428 gezelter 1126 for (bend = mol->beginBend(bendIter); bend != NULL;
429     bend = mol->nextBend(bendIter)) {
430    
431     RealType angle;
432     bend->calcForce(angle);
433     RealType currBendPot = bend->getPotential();
434 gezelter 1448
435 gezelter 1126 bendPotential += bend->getPotential();
436 gezelter 1545 map<Bend*, BendDataSet>::iterator i = bendDataSets.find(bend);
437 gezelter 1126 if (i == bendDataSets.end()) {
438     BendDataSet dataSet;
439     dataSet.prev.angle = dataSet.curr.angle = angle;
440     dataSet.prev.potential = dataSet.curr.potential = currBendPot;
441     dataSet.deltaV = 0.0;
442 gezelter 1590 bendDataSets.insert(map<Bend*, BendDataSet>::value_type(bend,
443     dataSet));
444 gezelter 1126 }else {
445     i->second.prev.angle = i->second.curr.angle;
446     i->second.prev.potential = i->second.curr.potential;
447     i->second.curr.angle = angle;
448     i->second.curr.potential = currBendPot;
449     i->second.deltaV = fabs(i->second.curr.potential -
450     i->second.prev.potential);
451     }
452 gezelter 507 }
453 gezelter 1126
454     for (torsion = mol->beginTorsion(torsionIter); torsion != NULL;
455     torsion = mol->nextTorsion(torsionIter)) {
456 tim 963 RealType angle;
457 gezelter 1126 torsion->calcForce(angle);
458 tim 963 RealType currTorsionPot = torsion->getPotential();
459 gezelter 1126 torsionPotential += torsion->getPotential();
460 gezelter 1545 map<Torsion*, TorsionDataSet>::iterator i = torsionDataSets.find(torsion);
461 gezelter 1126 if (i == torsionDataSets.end()) {
462     TorsionDataSet dataSet;
463     dataSet.prev.angle = dataSet.curr.angle = angle;
464     dataSet.prev.potential = dataSet.curr.potential = currTorsionPot;
465     dataSet.deltaV = 0.0;
466 gezelter 1545 torsionDataSets.insert(map<Torsion*, TorsionDataSet>::value_type(torsion, dataSet));
467 gezelter 1126 }else {
468     i->second.prev.angle = i->second.curr.angle;
469     i->second.prev.potential = i->second.curr.potential;
470     i->second.curr.angle = angle;
471     i->second.curr.potential = currTorsionPot;
472     i->second.deltaV = fabs(i->second.curr.potential -
473     i->second.prev.potential);
474     }
475     }
476 gezelter 1545
477 cli2 1275 for (inversion = mol->beginInversion(inversionIter);
478     inversion != NULL;
479     inversion = mol->nextInversion(inversionIter)) {
480     RealType angle;
481     inversion->calcForce(angle);
482     RealType currInversionPot = inversion->getPotential();
483     inversionPotential += inversion->getPotential();
484 gezelter 1545 map<Inversion*, InversionDataSet>::iterator i = inversionDataSets.find(inversion);
485 cli2 1275 if (i == inversionDataSets.end()) {
486     InversionDataSet dataSet;
487     dataSet.prev.angle = dataSet.curr.angle = angle;
488     dataSet.prev.potential = dataSet.curr.potential = currInversionPot;
489     dataSet.deltaV = 0.0;
490 gezelter 1545 inversionDataSets.insert(map<Inversion*, InversionDataSet>::value_type(inversion, dataSet));
491 cli2 1275 }else {
492     i->second.prev.angle = i->second.curr.angle;
493     i->second.prev.potential = i->second.curr.potential;
494     i->second.curr.angle = angle;
495     i->second.curr.potential = currInversionPot;
496     i->second.deltaV = fabs(i->second.curr.potential -
497     i->second.prev.potential);
498     }
499     }
500 gezelter 246 }
501    
502 gezelter 1126 RealType shortRangePotential = bondPotential + bendPotential +
503 cli2 1275 torsionPotential + inversionPotential;
504 gezelter 246 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
505     curSnapshot->statData[Stats::SHORT_RANGE_POTENTIAL] = shortRangePotential;
506 tim 665 curSnapshot->statData[Stats::BOND_POTENTIAL] = bondPotential;
507     curSnapshot->statData[Stats::BEND_POTENTIAL] = bendPotential;
508     curSnapshot->statData[Stats::DIHEDRAL_POTENTIAL] = torsionPotential;
509 gezelter 1545 curSnapshot->statData[Stats::INVERSION_POTENTIAL] = inversionPotential;
510 gezelter 507 }
511 gezelter 1126
512 chuckv 1595 void ForceManager::longRangeInteractionsRapaport() {
513    
514     Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
515     DataStorage* config = &(curSnapshot->atomData);
516     DataStorage* cgConfig = &(curSnapshot->cgData);
517    
518     //calculate the center of mass of cutoff group
519    
520     SimInfo::MoleculeIterator mi;
521     Molecule* mol;
522     Molecule::CutoffGroupIterator ci;
523     CutoffGroup* cg;
524    
525     if(info_->getNCutoffGroups() > 0){
526     for (mol = info_->beginMolecule(mi); mol != NULL;
527     mol = info_->nextMolecule(mi)) {
528     for(cg = mol->beginCutoffGroup(ci); cg != NULL;
529     cg = mol->nextCutoffGroup(ci)) {
530     cerr << "branch1\n";
531     cerr << "globind = " << cg->getGlobalIndex() << "\n";
532     cg->updateCOM();
533     }
534     }
535     } else {
536     // center of mass of the group is the same as position of the atom
537     // if cutoff group does not exist
538     cerr << "branch2\n";
539     cgConfig->position = config->position;
540     }
541    
542     fDecomp_->zeroWorkArrays();
543     fDecomp_->distributeData();
544    
545     int cg1, cg2, atom1, atom2, topoDist;
546     Vector3d d_grp, dag, d;
547     RealType rgrpsq, rgrp, r2, r;
548     RealType electroMult, vdwMult;
549     RealType vij;
550     Vector3d fij, fg, f1;
551     tuple3<RealType, RealType, RealType> cuts;
552     RealType rCutSq;
553     bool in_switching_region;
554     RealType sw, dswdr, swderiv;
555     vector<int> atomListColumn, atomListRow, atomListLocal;
556     InteractionData idat;
557     SelfData sdat;
558     RealType mf;
559     RealType lrPot;
560     RealType vpair;
561     potVec longRangePotential(0.0);
562     potVec workPot(0.0);
563    
564     int loopStart, loopEnd;
565    
566     idat.vdwMult = &vdwMult;
567     idat.electroMult = &electroMult;
568     idat.pot = &workPot;
569     sdat.pot = fDecomp_->getEmbeddingPotential();
570     idat.vpair = &vpair;
571     idat.f1 = &f1;
572     idat.sw = &sw;
573     idat.shiftedPot = (cutoffMethod_ == SHIFTED_POTENTIAL) ? true : false;
574     idat.shiftedForce = (cutoffMethod_ == SHIFTED_FORCE) ? true : false;
575    
576     loopEnd = PAIR_LOOP;
577     if (info_->requiresPrepair() ) {
578     loopStart = PREPAIR_LOOP;
579     } else {
580     loopStart = PAIR_LOOP;
581     }
582    
583     for (int iLoop = loopStart; iLoop <= loopEnd; iLoop++) {
584    
585     if (iLoop == loopStart) {
586     bool update_nlist = fDecomp_->checkNeighborList();
587     if (update_nlist)
588     neighborMatW = fDecomp_->buildLayerBasedNeighborList();
589     }
590    
591     int i;
592     #pragma omp parallel for num_threads(2) private(i)
593     for(i = 0; i < neighborMatW.size(); ++i)
594     for(vector<int>::iterator j = neighborMatW[i].begin(); j != neighborMatW[i].end(); ++j)
595     {
596     cg1 = i;
597     cg2 = *j;
598    
599     cuts = fDecomp_->getGroupCutoffs(cg1, cg2);
600    
601     d_grp = fDecomp_->getIntergroupVector(cg1, cg2);
602     curSnapshot->wrapVector(d_grp);
603     rgrpsq = d_grp.lengthSquare();
604    
605     rCutSq = cuts.second;
606    
607     if (rgrpsq < rCutSq) {
608     idat.rcut = &cuts.first;
609     if (iLoop == PAIR_LOOP) {
610     vij = 0.0;
611     fij = V3Zero;
612     }
613    
614     in_switching_region = switcher_->getSwitch(rgrpsq, sw, dswdr,
615     rgrp);
616    
617     atomListRow = fDecomp_->getAtomsInGroupRow(cg1);
618     atomListColumn = fDecomp_->getAtomsInGroupColumn(cg2);
619    
620     for (vector<int>::iterator ia = atomListRow.begin();
621     ia != atomListRow.end(); ++ia) {
622     atom1 = (*ia);
623    
624     for (vector<int>::iterator jb = atomListColumn.begin();
625     jb != atomListColumn.end(); ++jb) {
626     atom2 = (*jb);
627    
628     if (!fDecomp_->skipAtomPair(atom1, atom2)) {
629     vpair = 0.0;
630     workPot = 0.0;
631     f1 = V3Zero;
632    
633     fDecomp_->fillInteractionData(idat, atom1, atom2);
634    
635     topoDist = fDecomp_->getTopologicalDistance(atom1, atom2);
636     vdwMult = vdwScale_[topoDist];
637     electroMult = electrostaticScale_[topoDist];
638    
639     if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
640     idat.d = &d_grp;
641     idat.r2 = &rgrpsq;
642     cerr << "dgrp = " << d_grp << "\n";
643     } else {
644     d = fDecomp_->getInteratomicVector(atom1, atom2);
645     curSnapshot->wrapVector( d );
646     r2 = d.lengthSquare();
647     cerr << "datm = " << d<< "\n";
648     idat.d = &d;
649     idat.r2 = &r2;
650     }
651    
652     cerr << "idat.d = " << *(idat.d) << "\n";
653     r = sqrt( *(idat.r2) );
654     idat.rij = &r;
655    
656     if (iLoop == PREPAIR_LOOP) {
657     interactionMan_->doPrePair(idat);
658     } else {
659     interactionMan_->doPair(idat);
660     fDecomp_->unpackInteractionData(idat, atom1, atom2);
661    
662     cerr << "d = " << *(idat.d) << "\tv=" << vpair << "\tf=" << f1 << "\n";
663     vij += vpair;
664     fij += f1;
665     tau -= outProduct( *(idat.d), f1);
666     }
667     }
668     }
669     }
670    
671     if (iLoop == PAIR_LOOP) {
672     if (in_switching_region) {
673     swderiv = vij * dswdr / rgrp;
674     fg = swderiv * d_grp;
675     fij += fg;
676    
677     if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
678     tau -= outProduct( *(idat.d), fg);
679     }
680    
681     for (vector<int>::iterator ia = atomListRow.begin();
682     ia != atomListRow.end(); ++ia) {
683     atom1 = (*ia);
684     mf = fDecomp_->getMassFactorRow(atom1);
685     // fg is the force on atom ia due to cutoff group's
686     // presence in switching region
687     fg = swderiv * d_grp * mf;
688     fDecomp_->addForceToAtomRow(atom1, fg);
689    
690     if (atomListRow.size() > 1) {
691     if (info_->usesAtomicVirial()) {
692     // find the distance between the atom
693     // and the center of the cutoff group:
694     dag = fDecomp_->getAtomToGroupVectorRow(atom1, cg1);
695     tau -= outProduct(dag, fg);
696     }
697     }
698     }
699     for (vector<int>::iterator jb = atomListColumn.begin();
700     jb != atomListColumn.end(); ++jb) {
701     atom2 = (*jb);
702     mf = fDecomp_->getMassFactorColumn(atom2);
703     // fg is the force on atom jb due to cutoff group's
704     // presence in switching region
705     fg = -swderiv * d_grp * mf;
706     fDecomp_->addForceToAtomColumn(atom2, fg);
707    
708     if (atomListColumn.size() > 1) {
709     if (info_->usesAtomicVirial()) {
710     // find the distance between the atom
711     // and the center of the cutoff group:
712     dag = fDecomp_->getAtomToGroupVectorColumn(atom2, cg2);
713     tau -= outProduct(dag, fg);
714     }
715     }
716     }
717     }
718     //if (!SIM_uses_AtomicVirial) {
719     // tau -= outProduct(d_grp, fij);
720     //}
721     }
722     }
723     }
724    
725     if (iLoop == PREPAIR_LOOP) {
726     if (info_->requiresPrepair()) {
727    
728     fDecomp_->collectIntermediateData();
729    
730     for (int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
731     fDecomp_->fillSelfData(sdat, atom1);
732     interactionMan_->doPreForce(sdat);
733     }
734    
735     fDecomp_->distributeIntermediateData();
736    
737     }
738     }
739    
740     }
741    
742     fDecomp_->collectData();
743    
744     if (info_->requiresSelfCorrection()) {
745    
746     for (int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
747     fDecomp_->fillSelfData(sdat, atom1);
748     interactionMan_->doSelfCorrection(sdat);
749     }
750    
751     }
752    
753     longRangePotential = *(fDecomp_->getEmbeddingPotential()) +
754     *(fDecomp_->getPairwisePotential());
755    
756     lrPot = longRangePotential.sum();
757    
758     //store the tau and long range potential
759     curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL] = lrPot;
760     curSnapshot->statData[Stats::VANDERWAALS_POTENTIAL] = longRangePotential[VANDERWAALS_FAMILY];
761     curSnapshot->statData[Stats::ELECTROSTATIC_POTENTIAL] = longRangePotential[ELECTROSTATIC_FAMILY];
762     }
763    
764 gezelter 1546 void ForceManager::longRangeInteractions() {
765 gezelter 1581
766 gezelter 1545 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
767     DataStorage* config = &(curSnapshot->atomData);
768     DataStorage* cgConfig = &(curSnapshot->cgData);
769    
770 gezelter 1581 //calculate the center of mass of cutoff group
771    
772     SimInfo::MoleculeIterator mi;
773     Molecule* mol;
774     Molecule::CutoffGroupIterator ci;
775     CutoffGroup* cg;
776    
777     if(info_->getNCutoffGroups() > 0){
778     for (mol = info_->beginMolecule(mi); mol != NULL;
779     mol = info_->nextMolecule(mi)) {
780     for(cg = mol->beginCutoffGroup(ci); cg != NULL;
781     cg = mol->nextCutoffGroup(ci)) {
782 gezelter 1593 cerr << "branch1\n";
783     cerr << "globind = " << cg->getGlobalIndex() << "\n";
784 gezelter 1581 cg->updateCOM();
785     }
786     }
787     } else {
788     // center of mass of the group is the same as position of the atom
789     // if cutoff group does not exist
790 gezelter 1593 cerr << "branch2\n";
791 gezelter 1581 cgConfig->position = config->position;
792     }
793    
794 gezelter 1575 fDecomp_->zeroWorkArrays();
795 gezelter 1549 fDecomp_->distributeData();
796 gezelter 1579
797     int cg1, cg2, atom1, atom2, topoDist;
798     Vector3d d_grp, dag, d;
799     RealType rgrpsq, rgrp, r2, r;
800     RealType electroMult, vdwMult;
801 gezelter 1549 RealType vij;
802 gezelter 1581 Vector3d fij, fg, f1;
803 gezelter 1576 tuple3<RealType, RealType, RealType> cuts;
804 gezelter 1545 RealType rCutSq;
805     bool in_switching_region;
806     RealType sw, dswdr, swderiv;
807 gezelter 1549 vector<int> atomListColumn, atomListRow, atomListLocal;
808 gezelter 1545 InteractionData idat;
809 gezelter 1546 SelfData sdat;
810     RealType mf;
811 gezelter 1575 RealType lrPot;
812 gezelter 1579 RealType vpair;
813 gezelter 1583 potVec longRangePotential(0.0);
814     potVec workPot(0.0);
815 gezelter 1544
816 gezelter 1545 int loopStart, loopEnd;
817 gezelter 1544
818 gezelter 1581 idat.vdwMult = &vdwMult;
819     idat.electroMult = &electroMult;
820 gezelter 1583 idat.pot = &workPot;
821     sdat.pot = fDecomp_->getEmbeddingPotential();
822 gezelter 1581 idat.vpair = &vpair;
823     idat.f1 = &f1;
824     idat.sw = &sw;
825 gezelter 1583 idat.shiftedPot = (cutoffMethod_ == SHIFTED_POTENTIAL) ? true : false;
826     idat.shiftedForce = (cutoffMethod_ == SHIFTED_FORCE) ? true : false;
827    
828 gezelter 1545 loopEnd = PAIR_LOOP;
829 gezelter 1546 if (info_->requiresPrepair() ) {
830 gezelter 1545 loopStart = PREPAIR_LOOP;
831     } else {
832     loopStart = PAIR_LOOP;
833     }
834 gezelter 1583
835 gezelter 1579 for (int iLoop = loopStart; iLoop <= loopEnd; iLoop++) {
836    
837 gezelter 1545 if (iLoop == loopStart) {
838 gezelter 1549 bool update_nlist = fDecomp_->checkNeighborList();
839 chuckv 1595 if (update_nlist)
840 gezelter 1549 neighborList = fDecomp_->buildNeighborList();
841 chuckv 1595
842 gezelter 1579 }
843    
844 chuckv 1595 for (vector<pair<int, int> >::iterator it = neighborList.begin();
845     it != neighborList.end(); ++it)
846     {
847 gezelter 1545 cg1 = (*it).first;
848     cg2 = (*it).second;
849 gezelter 1576
850     cuts = fDecomp_->getGroupCutoffs(cg1, cg2);
851 gezelter 1545
852 gezelter 1549 d_grp = fDecomp_->getIntergroupVector(cg1, cg2);
853 gezelter 1545 curSnapshot->wrapVector(d_grp);
854     rgrpsq = d_grp.lengthSquare();
855    
856 gezelter 1576 rCutSq = cuts.second;
857    
858 gezelter 1545 if (rgrpsq < rCutSq) {
859 gezelter 1579 idat.rcut = &cuts.first;
860 gezelter 1545 if (iLoop == PAIR_LOOP) {
861 gezelter 1587 vij = 0.0;
862 gezelter 1545 fij = V3Zero;
863     }
864    
865 gezelter 1579 in_switching_region = switcher_->getSwitch(rgrpsq, sw, dswdr,
866 gezelter 1576 rgrp);
867    
868 gezelter 1549 atomListRow = fDecomp_->getAtomsInGroupRow(cg1);
869     atomListColumn = fDecomp_->getAtomsInGroupColumn(cg2);
870 gezelter 1545
871 gezelter 1549 for (vector<int>::iterator ia = atomListRow.begin();
872     ia != atomListRow.end(); ++ia) {
873 gezelter 1545 atom1 = (*ia);
874    
875 gezelter 1549 for (vector<int>::iterator jb = atomListColumn.begin();
876     jb != atomListColumn.end(); ++jb) {
877 gezelter 1545 atom2 = (*jb);
878 gezelter 1593
879 gezelter 1549 if (!fDecomp_->skipAtomPair(atom1, atom2)) {
880 gezelter 1579 vpair = 0.0;
881 gezelter 1583 workPot = 0.0;
882 gezelter 1581 f1 = V3Zero;
883 gezelter 1575
884 gezelter 1581 fDecomp_->fillInteractionData(idat, atom1, atom2);
885 gezelter 1579
886     topoDist = fDecomp_->getTopologicalDistance(atom1, atom2);
887     vdwMult = vdwScale_[topoDist];
888     electroMult = electrostaticScale_[topoDist];
889 gezelter 1546
890 gezelter 1549 if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
891 gezelter 1579 idat.d = &d_grp;
892     idat.r2 = &rgrpsq;
893 gezelter 1593 cerr << "dgrp = " << d_grp << "\n";
894 gezelter 1545 } else {
895 gezelter 1579 d = fDecomp_->getInteratomicVector(atom1, atom2);
896     curSnapshot->wrapVector( d );
897     r2 = d.lengthSquare();
898 gezelter 1593 cerr << "datm = " << d<< "\n";
899 gezelter 1579 idat.d = &d;
900     idat.r2 = &r2;
901 gezelter 1545 }
902    
903 gezelter 1593 cerr << "idat.d = " << *(idat.d) << "\n";
904 gezelter 1581 r = sqrt( *(idat.r2) );
905 gezelter 1579 idat.rij = &r;
906 gezelter 1546
907 gezelter 1545 if (iLoop == PREPAIR_LOOP) {
908     interactionMan_->doPrePair(idat);
909     } else {
910     interactionMan_->doPair(idat);
911 gezelter 1575 fDecomp_->unpackInteractionData(idat, atom1, atom2);
912 gezelter 1593
913     cerr << "d = " << *(idat.d) << "\tv=" << vpair << "\tf=" << f1 << "\n";
914 gezelter 1581 vij += vpair;
915     fij += f1;
916     tau -= outProduct( *(idat.d), f1);
917 gezelter 1545 }
918     }
919     }
920     }
921    
922     if (iLoop == PAIR_LOOP) {
923     if (in_switching_region) {
924     swderiv = vij * dswdr / rgrp;
925     fg = swderiv * d_grp;
926     fij += fg;
927    
928 gezelter 1549 if (atomListRow.size() == 1 && atomListColumn.size() == 1) {
929 gezelter 1554 tau -= outProduct( *(idat.d), fg);
930 gezelter 1545 }
931    
932 gezelter 1549 for (vector<int>::iterator ia = atomListRow.begin();
933     ia != atomListRow.end(); ++ia) {
934 gezelter 1545 atom1 = (*ia);
935 gezelter 1569 mf = fDecomp_->getMassFactorRow(atom1);
936 gezelter 1545 // fg is the force on atom ia due to cutoff group's
937     // presence in switching region
938     fg = swderiv * d_grp * mf;
939 gezelter 1549 fDecomp_->addForceToAtomRow(atom1, fg);
940 gezelter 1545
941 gezelter 1549 if (atomListRow.size() > 1) {
942 gezelter 1546 if (info_->usesAtomicVirial()) {
943 gezelter 1545 // find the distance between the atom
944     // and the center of the cutoff group:
945 gezelter 1549 dag = fDecomp_->getAtomToGroupVectorRow(atom1, cg1);
946 gezelter 1545 tau -= outProduct(dag, fg);
947     }
948     }
949     }
950 gezelter 1549 for (vector<int>::iterator jb = atomListColumn.begin();
951     jb != atomListColumn.end(); ++jb) {
952 gezelter 1545 atom2 = (*jb);
953 gezelter 1569 mf = fDecomp_->getMassFactorColumn(atom2);
954 gezelter 1545 // fg is the force on atom jb due to cutoff group's
955     // presence in switching region
956     fg = -swderiv * d_grp * mf;
957 gezelter 1549 fDecomp_->addForceToAtomColumn(atom2, fg);
958 gezelter 1545
959 gezelter 1549 if (atomListColumn.size() > 1) {
960 gezelter 1546 if (info_->usesAtomicVirial()) {
961 gezelter 1545 // find the distance between the atom
962     // and the center of the cutoff group:
963 gezelter 1549 dag = fDecomp_->getAtomToGroupVectorColumn(atom2, cg2);
964 gezelter 1545 tau -= outProduct(dag, fg);
965     }
966     }
967     }
968     }
969     //if (!SIM_uses_AtomicVirial) {
970     // tau -= outProduct(d_grp, fij);
971     //}
972     }
973     }
974     }
975    
976     if (iLoop == PREPAIR_LOOP) {
977 gezelter 1590 if (info_->requiresPrepair()) {
978    
979 gezelter 1549 fDecomp_->collectIntermediateData();
980 gezelter 1570
981     for (int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
982 gezelter 1581 fDecomp_->fillSelfData(sdat, atom1);
983 gezelter 1545 interactionMan_->doPreForce(sdat);
984     }
985 gezelter 1590
986     fDecomp_->distributeIntermediateData();
987    
988 gezelter 1545 }
989     }
990    
991 gezelter 1544 }
992 gezelter 1545
993 gezelter 1549 fDecomp_->collectData();
994 gezelter 1570
995     if (info_->requiresSelfCorrection()) {
996 gezelter 1545
997 gezelter 1570 for (int atom1 = 0; atom1 < info_->getNAtoms(); atom1++) {
998 gezelter 1581 fDecomp_->fillSelfData(sdat, atom1);
999 gezelter 1570 interactionMan_->doSelfCorrection(sdat);
1000     }
1001    
1002     }
1003    
1004 gezelter 1583 longRangePotential = *(fDecomp_->getEmbeddingPotential()) +
1005     *(fDecomp_->getPairwisePotential());
1006    
1007 gezelter 1575 lrPot = longRangePotential.sum();
1008    
1009 gezelter 246 //store the tau and long range potential
1010 chuckv 664 curSnapshot->statData[Stats::LONG_RANGE_POTENTIAL] = lrPot;
1011 gezelter 1550 curSnapshot->statData[Stats::VANDERWAALS_POTENTIAL] = longRangePotential[VANDERWAALS_FAMILY];
1012     curSnapshot->statData[Stats::ELECTROSTATIC_POTENTIAL] = longRangePotential[ELECTROSTATIC_FAMILY];
1013 gezelter 507 }
1014 gezelter 246
1015 gezelter 1126
1016 gezelter 1464 void ForceManager::postCalculation() {
1017 gezelter 246 SimInfo::MoleculeIterator mi;
1018     Molecule* mol;
1019     Molecule::RigidBodyIterator rbIter;
1020     RigidBody* rb;
1021 gezelter 1126 Snapshot* curSnapshot = info_->getSnapshotManager()->getCurrentSnapshot();
1022 gezelter 246
1023     // collect the atomic forces onto rigid bodies
1024 gezelter 1126
1025     for (mol = info_->beginMolecule(mi); mol != NULL;
1026     mol = info_->nextMolecule(mi)) {
1027     for (rb = mol->beginRigidBody(rbIter); rb != NULL;
1028     rb = mol->nextRigidBody(rbIter)) {
1029 gezelter 1464 Mat3x3d rbTau = rb->calcForcesAndTorquesAndVirial();
1030     tau += rbTau;
1031 gezelter 507 }
1032 gezelter 1126 }
1033 gezelter 1464
1034 gezelter 1126 #ifdef IS_MPI
1035 gezelter 1464 Mat3x3d tmpTau(tau);
1036     MPI_Allreduce(tmpTau.getArrayPointer(), tau.getArrayPointer(),
1037     9, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
1038 gezelter 1126 #endif
1039 gezelter 1464 curSnapshot->statData.setTau(tau);
1040 gezelter 507 }
1041 gezelter 246
1042 gezelter 1390 } //end namespace OpenMD

Properties

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