ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/ForceManager.cpp
Revision: 2060
Committed: Tue Mar 3 16:19:47 2015 UTC (10 years, 1 month ago) by gezelter
File size: 36968 byte(s)
Log Message:
Got rid of unused variables and spurious warnings

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

Properties

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