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