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