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] Vardeman & Gezelter, in progress (2009). |
40 |
*/ |
41 |
|
42 |
#include "nonbonded/InteractionManager.hpp" |
43 |
#include "UseTheForce/doForces_interface.h" |
44 |
|
45 |
namespace OpenMD { |
46 |
|
47 |
InteractionManager* InteractionManager::_instance = NULL; |
48 |
SimInfo* InteractionManager::info_ = NULL; |
49 |
bool InteractionManager::initialized_ = false; |
50 |
|
51 |
RealType InteractionManager::rCut_ = 0.0; |
52 |
RealType InteractionManager::rSwitch_ = 0.0; |
53 |
RealType InteractionManager::skinThickness_ = 0.0; |
54 |
RealType InteractionManager::listRadius_ = 0.0; |
55 |
CutoffMethod InteractionManager::cutoffMethod_ = SHIFTED_FORCE; |
56 |
SwitchingFunctionType InteractionManager::sft_ = cubic; |
57 |
RealType InteractionManager::vdwScale_[4] = {1.0, 0.0, 0.0, 0.0}; |
58 |
RealType InteractionManager::electrostaticScale_[4] = {1.0, 0.0, 0.0, 0.0}; |
59 |
|
60 |
map<int, AtomType*> InteractionManager::typeMap_; |
61 |
map<pair<AtomType*, AtomType*>, set<NonBondedInteraction*> > InteractionManager::interactions_; |
62 |
|
63 |
LJ* InteractionManager::lj_ = new LJ(); |
64 |
GB* InteractionManager::gb_ = new GB(); |
65 |
Sticky* InteractionManager::sticky_ = new Sticky(); |
66 |
Morse* InteractionManager::morse_ = new Morse(); |
67 |
EAM* InteractionManager::eam_ = new EAM(); |
68 |
SC* InteractionManager::sc_ = new SC(); |
69 |
Electrostatic* InteractionManager::electrostatic_ = new Electrostatic(); |
70 |
MAW* InteractionManager::maw_ = new MAW(); |
71 |
SwitchingFunction* InteractionManager::switcher_ = new SwitchingFunction(); |
72 |
|
73 |
InteractionManager* InteractionManager::Instance() { |
74 |
if (!_instance) { |
75 |
_instance = new InteractionManager(); |
76 |
} |
77 |
return _instance; |
78 |
} |
79 |
|
80 |
void InteractionManager::initialize() { |
81 |
|
82 |
ForceField* forceField_ = info_->getForceField(); |
83 |
|
84 |
lj_->setForceField(forceField_); |
85 |
gb_->setForceField(forceField_); |
86 |
sticky_->setForceField(forceField_); |
87 |
eam_->setForceField(forceField_); |
88 |
sc_->setForceField(forceField_); |
89 |
morse_->setForceField(forceField_); |
90 |
electrostatic_->setForceField(forceField_); |
91 |
maw_->setForceField(forceField_); |
92 |
|
93 |
ForceFieldOptions& fopts = forceField_->getForceFieldOptions(); |
94 |
|
95 |
// Force fields can set options on how to scale van der Waals and electrostatic |
96 |
// interactions for atoms connected via bonds, bends and torsions |
97 |
// in this case the topological distance between atoms is: |
98 |
// 0 = topologically unconnected |
99 |
// 1 = bonded together |
100 |
// 2 = connected via a bend |
101 |
// 3 = connected via a torsion |
102 |
|
103 |
vdwScale_[0] = 1.0; |
104 |
vdwScale_[1] = fopts.getvdw12scale(); |
105 |
vdwScale_[2] = fopts.getvdw13scale(); |
106 |
vdwScale_[3] = fopts.getvdw14scale(); |
107 |
|
108 |
electrostaticScale_[0] = 1.0; |
109 |
electrostaticScale_[1] = fopts.getelectrostatic12scale(); |
110 |
electrostaticScale_[2] = fopts.getelectrostatic13scale(); |
111 |
electrostaticScale_[3] = fopts.getelectrostatic14scale(); |
112 |
|
113 |
ForceField::AtomTypeContainer* atomTypes = forceField_->getAtomTypes(); |
114 |
ForceField::AtomTypeContainer::MapTypeIterator i1, i2; |
115 |
AtomType* atype1; |
116 |
AtomType* atype2; |
117 |
pair<AtomType*, AtomType*> key; |
118 |
pair<set<NonBondedInteraction*>::iterator, bool> ret; |
119 |
|
120 |
for (atype1 = atomTypes->beginType(i1); atype1 != NULL; |
121 |
atype1 = atomTypes->nextType(i1)) { |
122 |
|
123 |
// add it to the map: |
124 |
AtomTypeProperties atp = atype1->getATP(); |
125 |
|
126 |
pair<map<int,AtomType*>::iterator,bool> ret; |
127 |
ret = typeMap_.insert( pair<int, AtomType*>(atp.ident, atype1) ); |
128 |
if (ret.second == false) { |
129 |
sprintf( painCave.errMsg, |
130 |
"InteractionManager already had a previous entry with ident %d\n", |
131 |
atp.ident); |
132 |
painCave.severity = OPENMD_INFO; |
133 |
painCave.isFatal = 0; |
134 |
simError(); |
135 |
} |
136 |
} |
137 |
|
138 |
// Now, iterate over all known types and add to the interaction map: |
139 |
|
140 |
map<int, AtomType*>::iterator it1, it2; |
141 |
for (it1 = typeMap_.begin(); it1 != typeMap_.end(); ++it1) { |
142 |
atype1 = (*it1).second; |
143 |
|
144 |
for( it2 = typeMap_.begin(); it2 != typeMap_.end(); ++it2) { |
145 |
atype2 = (*it2).second; |
146 |
|
147 |
bool vdwExplicit = false; |
148 |
bool metExplicit = false; |
149 |
bool hbExplicit = false; |
150 |
|
151 |
key = make_pair(atype1, atype2); |
152 |
|
153 |
if (atype1->isLennardJones() && atype2->isLennardJones()) { |
154 |
interactions_[key].insert(lj_); |
155 |
} |
156 |
if (atype1->isElectrostatic() && atype2->isElectrostatic() ) { |
157 |
interactions_[key].insert(electrostatic_); |
158 |
} |
159 |
if (atype1->isSticky() && atype2->isSticky() ) { |
160 |
interactions_[key].insert(sticky_); |
161 |
} |
162 |
if (atype1->isStickyPower() && atype2->isStickyPower() ) { |
163 |
interactions_[key].insert(sticky_); |
164 |
} |
165 |
if (atype1->isEAM() && atype2->isEAM() ) { |
166 |
interactions_[key].insert(eam_); |
167 |
} |
168 |
if (atype1->isSC() && atype2->isSC() ) { |
169 |
interactions_[key].insert(sc_); |
170 |
} |
171 |
if (atype1->isGayBerne() && atype2->isGayBerne() ) { |
172 |
interactions_[key].insert(gb_); |
173 |
} |
174 |
if ((atype1->isGayBerne() && atype2->isLennardJones()) |
175 |
|| (atype1->isLennardJones() && atype2->isGayBerne())) { |
176 |
interactions_[key].insert(gb_); |
177 |
} |
178 |
|
179 |
// look for an explicitly-set non-bonded interaction type using the |
180 |
// two atom types. |
181 |
NonBondedInteractionType* nbiType = forceField_->getNonBondedInteractionType(atype1->getName(), atype2->getName()); |
182 |
|
183 |
if (nbiType != NULL) { |
184 |
|
185 |
if (nbiType->isLennardJones()) { |
186 |
// We found an explicit Lennard-Jones interaction. |
187 |
// override all other vdw entries for this pair of atom types: |
188 |
set<NonBondedInteraction*>::iterator it; |
189 |
for (it = interactions_[key].begin(); |
190 |
it != interactions_[key].end(); ++it) { |
191 |
InteractionFamily ifam = (*it)->getFamily(); |
192 |
if (ifam == VANDERWAALS_FAMILY) interactions_[key].erase(*it); |
193 |
} |
194 |
interactions_[key].insert(lj_); |
195 |
vdwExplicit = true; |
196 |
} |
197 |
|
198 |
if (nbiType->isMorse()) { |
199 |
if (vdwExplicit) { |
200 |
sprintf( painCave.errMsg, |
201 |
"InteractionManager::initialize found more than one " |
202 |
"explicit \n" |
203 |
"\tvan der Waals interaction for atom types %s - %s\n", |
204 |
atype1->getName().c_str(), atype2->getName().c_str()); |
205 |
painCave.severity = OPENMD_ERROR; |
206 |
painCave.isFatal = 1; |
207 |
simError(); |
208 |
} |
209 |
// We found an explicit Morse interaction. |
210 |
// override all other vdw entries for this pair of atom types: |
211 |
set<NonBondedInteraction*>::iterator it; |
212 |
for (it = interactions_[key].begin(); |
213 |
it != interactions_[key].end(); ++it) { |
214 |
InteractionFamily ifam = (*it)->getFamily(); |
215 |
if (ifam == VANDERWAALS_FAMILY) interactions_[key].erase(*it); |
216 |
} |
217 |
interactions_[key].insert(morse_); |
218 |
vdwExplicit = true; |
219 |
} |
220 |
|
221 |
if (nbiType->isEAM()) { |
222 |
// We found an explicit EAM interaction. |
223 |
// override all other metallic entries for this pair of atom types: |
224 |
set<NonBondedInteraction*>::iterator it; |
225 |
for (it = interactions_[key].begin(); |
226 |
it != interactions_[key].end(); ++it) { |
227 |
InteractionFamily ifam = (*it)->getFamily(); |
228 |
if (ifam == METALLIC_FAMILY) interactions_[key].erase(*it); |
229 |
} |
230 |
interactions_[key].insert(eam_); |
231 |
metExplicit = true; |
232 |
} |
233 |
|
234 |
if (nbiType->isSC()) { |
235 |
if (metExplicit) { |
236 |
sprintf( painCave.errMsg, |
237 |
"InteractionManager::initialize found more than one " |
238 |
"explicit\n" |
239 |
"\tmetallic interaction for atom types %s - %s\n", |
240 |
atype1->getName().c_str(), atype2->getName().c_str()); |
241 |
painCave.severity = OPENMD_ERROR; |
242 |
painCave.isFatal = 1; |
243 |
simError(); |
244 |
} |
245 |
// We found an explicit Sutton-Chen interaction. |
246 |
// override all other metallic entries for this pair of atom types: |
247 |
set<NonBondedInteraction*>::iterator it; |
248 |
for (it = interactions_[key].begin(); |
249 |
it != interactions_[key].end(); ++it) { |
250 |
InteractionFamily ifam = (*it)->getFamily(); |
251 |
if (ifam == METALLIC_FAMILY) interactions_[key].erase(*it); |
252 |
} |
253 |
interactions_[key].insert(sc_); |
254 |
metExplicit = true; |
255 |
} |
256 |
|
257 |
if (nbiType->isMAW()) { |
258 |
if (vdwExplicit) { |
259 |
sprintf( painCave.errMsg, |
260 |
"InteractionManager::initialize found more than one " |
261 |
"explicit\n" |
262 |
"\tvan der Waals interaction for atom types %s - %s\n", |
263 |
atype1->getName().c_str(), atype2->getName().c_str()); |
264 |
painCave.severity = OPENMD_ERROR; |
265 |
painCave.isFatal = 1; |
266 |
simError(); |
267 |
} |
268 |
// We found an explicit MAW interaction. |
269 |
// override all other vdw entries for this pair of atom types: |
270 |
set<NonBondedInteraction*>::iterator it; |
271 |
for (it = interactions_[key].begin(); |
272 |
it != interactions_[key].end(); ++it) { |
273 |
InteractionFamily ifam = (*it)->getFamily(); |
274 |
if (ifam == VANDERWAALS_FAMILY) interactions_[key].erase(*it); |
275 |
} |
276 |
interactions_[key].insert(maw_); |
277 |
vdwExplicit = true; |
278 |
} |
279 |
} |
280 |
} |
281 |
} |
282 |
|
283 |
|
284 |
// make sure every pair of atom types in this simulation has a |
285 |
// non-bonded interaction: |
286 |
|
287 |
set<AtomType*> simTypes = info_->getSimulatedAtomTypes(); |
288 |
set<AtomType*>::iterator it, jt; |
289 |
for (it = simTypes.begin(); it != simTypes.end(); ++it) { |
290 |
atype1 = (*it); |
291 |
for (jt = simTypes.begin(); jt != simTypes.end(); ++jt) { |
292 |
atype2 = (*jt); |
293 |
key = make_pair(atype1, atype2); |
294 |
|
295 |
if (interactions_[key].size() == 0) { |
296 |
sprintf( painCave.errMsg, |
297 |
"InteractionManager unable to find an appropriate non-bonded\n" |
298 |
"\tinteraction for atom types %s - %s\n", |
299 |
atype1->getName().c_str(), atype2->getName().c_str()); |
300 |
painCave.severity = OPENMD_INFO; |
301 |
painCave.isFatal = 1; |
302 |
simError(); |
303 |
} |
304 |
} |
305 |
} |
306 |
|
307 |
setupCutoffs(); |
308 |
setupSwitching(); |
309 |
setupNeighborlists(); |
310 |
|
311 |
int ljsp = cutoffMethod_ == SHIFTED_POTENTIAL ? 1 : 0; |
312 |
int ljsf = cutoffMethod_ == SHIFTED_FORCE ? 1 : 0; |
313 |
notifyFortranCutoffs(&rCut_, &rSwitch_, &ljsp, &ljsf); |
314 |
notifyFortranSkinThickness(&skinThickness_); |
315 |
|
316 |
initialized_ = true; |
317 |
} |
318 |
|
319 |
/** |
320 |
* setupCutoffs |
321 |
* |
322 |
* Sets the values of cutoffRadius and cutoffMethod |
323 |
* |
324 |
* cutoffRadius : realType |
325 |
* If the cutoffRadius was explicitly set, use that value. |
326 |
* If the cutoffRadius was not explicitly set: |
327 |
* Are there electrostatic atoms? Use 12.0 Angstroms. |
328 |
* No electrostatic atoms? Poll the atom types present in the |
329 |
* simulation for suggested cutoff values (e.g. 2.5 * sigma). |
330 |
* Use the maximum suggested value that was found. |
331 |
* |
332 |
* cutoffMethod : (one of HARD, SWITCHED, SHIFTED_FORCE, SHIFTED_POTENTIAL) |
333 |
* If cutoffMethod was explicitly set, use that choice. |
334 |
* If cutoffMethod was not explicitly set, use SHIFTED_FORCE |
335 |
*/ |
336 |
void InteractionManager::setupCutoffs() { |
337 |
|
338 |
Globals* simParams_ = info_->getSimParams(); |
339 |
|
340 |
if (simParams_->haveCutoffRadius()) { |
341 |
rCut_ = simParams_->getCutoffRadius(); |
342 |
} else { |
343 |
if (info_->usesElectrostaticAtoms()) { |
344 |
sprintf(painCave.errMsg, |
345 |
"InteractionManager::setupCutoffs: No value was set for the cutoffRadius.\n" |
346 |
"\tOpenMD will use a default value of 12.0 angstroms" |
347 |
"\tfor the cutoffRadius.\n"); |
348 |
painCave.isFatal = 0; |
349 |
painCave.severity = OPENMD_INFO; |
350 |
simError(); |
351 |
rCut_ = 12.0; |
352 |
} else { |
353 |
RealType thisCut; |
354 |
set<AtomType*>::iterator i; |
355 |
set<AtomType*> atomTypes; |
356 |
atomTypes = info_->getSimulatedAtomTypes(); |
357 |
for (i = atomTypes.begin(); i != atomTypes.end(); ++i) { |
358 |
thisCut = getSuggestedCutoffRadius((*i)); |
359 |
rCut_ = max(thisCut, rCut_); |
360 |
} |
361 |
sprintf(painCave.errMsg, |
362 |
"InteractionManager::setupCutoffs: No value was set for the cutoffRadius.\n" |
363 |
"\tOpenMD will use %lf angstroms.\n", |
364 |
rCut_); |
365 |
painCave.isFatal = 0; |
366 |
painCave.severity = OPENMD_INFO; |
367 |
simError(); |
368 |
} |
369 |
} |
370 |
|
371 |
map<string, CutoffMethod> stringToCutoffMethod; |
372 |
stringToCutoffMethod["HARD"] = HARD; |
373 |
stringToCutoffMethod["SWITCHED"] = SWITCHED; |
374 |
stringToCutoffMethod["SHIFTED_POTENTIAL"] = SHIFTED_POTENTIAL; |
375 |
stringToCutoffMethod["SHIFTED_FORCE"] = SHIFTED_FORCE; |
376 |
|
377 |
if (simParams_->haveCutoffMethod()) { |
378 |
string cutMeth = toUpperCopy(simParams_->getCutoffMethod()); |
379 |
map<string, CutoffMethod>::iterator i; |
380 |
i = stringToCutoffMethod.find(cutMeth); |
381 |
if (i == stringToCutoffMethod.end()) { |
382 |
sprintf(painCave.errMsg, |
383 |
"InteractionManager::setupCutoffs: Could not find chosen cutoffMethod %s\n" |
384 |
"\tShould be one of: " |
385 |
"HARD, SWITCHED, SHIFTED_POTENTIAL, or SHIFTED_FORCE\n", |
386 |
cutMeth.c_str()); |
387 |
painCave.isFatal = 1; |
388 |
painCave.severity = OPENMD_ERROR; |
389 |
simError(); |
390 |
} else { |
391 |
cutoffMethod_ = i->second; |
392 |
} |
393 |
} else { |
394 |
sprintf(painCave.errMsg, |
395 |
"InteractionManager::setupCutoffs: No value was set for the cutoffMethod.\n" |
396 |
"\tOpenMD will use SHIFTED_FORCE.\n"); |
397 |
painCave.isFatal = 0; |
398 |
painCave.severity = OPENMD_INFO; |
399 |
simError(); |
400 |
cutoffMethod_ = SHIFTED_FORCE; |
401 |
} |
402 |
} |
403 |
|
404 |
|
405 |
/** |
406 |
* setupSwitching |
407 |
* |
408 |
* Sets the values of switchingRadius and |
409 |
* If the switchingRadius was explicitly set, use that value (but check it) |
410 |
* If the switchingRadius was not explicitly set: use 0.85 * cutoffRadius_ |
411 |
*/ |
412 |
void InteractionManager::setupSwitching() { |
413 |
Globals* simParams_ = info_->getSimParams(); |
414 |
|
415 |
if (simParams_->haveSwitchingRadius()) { |
416 |
rSwitch_ = simParams_->getSwitchingRadius(); |
417 |
if (rSwitch_ > rCut_) { |
418 |
sprintf(painCave.errMsg, |
419 |
"InteractionManager::setupSwitching: switchingRadius (%f) is larger than cutoffRadius(%f)\n", |
420 |
rSwitch_, rCut_); |
421 |
painCave.isFatal = 1; |
422 |
painCave.severity = OPENMD_ERROR; |
423 |
simError(); |
424 |
} |
425 |
} else { |
426 |
rSwitch_ = 0.85 * rCut_; |
427 |
sprintf(painCave.errMsg, |
428 |
"InteractionManager::setupSwitching: No value was set for the switchingRadius.\n" |
429 |
"\tOpenMD will use a default value of 85 percent of the cutoffRadius.\n" |
430 |
"\tswitchingRadius = %f. for this simulation\n", rSwitch_); |
431 |
painCave.isFatal = 0; |
432 |
painCave.severity = OPENMD_WARNING; |
433 |
simError(); |
434 |
} |
435 |
|
436 |
if (simParams_->haveSwitchingFunctionType()) { |
437 |
string funcType = simParams_->getSwitchingFunctionType(); |
438 |
toUpper(funcType); |
439 |
if (funcType == "CUBIC") { |
440 |
sft_ = cubic; |
441 |
} else { |
442 |
if (funcType == "FIFTH_ORDER_POLYNOMIAL") { |
443 |
sft_ = fifth_order_poly; |
444 |
} else { |
445 |
// throw error |
446 |
sprintf( painCave.errMsg, |
447 |
"InteractionManager::setupSwitching : Unknown switchingFunctionType. (Input file specified %s .)\n" |
448 |
"\tswitchingFunctionType must be one of: " |
449 |
"\"cubic\" or \"fifth_order_polynomial\".", |
450 |
funcType.c_str() ); |
451 |
painCave.isFatal = 1; |
452 |
painCave.severity = OPENMD_ERROR; |
453 |
simError(); |
454 |
} |
455 |
} |
456 |
} |
457 |
|
458 |
switcher_->setSwitchType(sft_); |
459 |
switcher_->setSwitch(rSwitch_, rCut_); |
460 |
} |
461 |
|
462 |
/** |
463 |
* setupNeighborlists |
464 |
* |
465 |
* If the skinThickness was explicitly set, use that value (but check it) |
466 |
* If the skinThickness was not explicitly set: use 1.0 angstroms |
467 |
*/ |
468 |
void InteractionManager::setupNeighborlists() { |
469 |
|
470 |
Globals* simParams_ = info_->getSimParams(); |
471 |
|
472 |
if (simParams_->haveSkinThickness()) { |
473 |
skinThickness_ = simParams_->getSkinThickness(); |
474 |
} else { |
475 |
skinThickness_ = 1.0; |
476 |
sprintf(painCave.errMsg, |
477 |
"InteractionManager::setupNeighborlists: No value was set for the skinThickness.\n" |
478 |
"\tOpenMD will use a default value of %f Angstroms\n" |
479 |
"\tfor this simulation\n", skinThickness_); |
480 |
painCave.severity = OPENMD_INFO; |
481 |
painCave.isFatal = 0; |
482 |
simError(); |
483 |
} |
484 |
|
485 |
listRadius_ = rCut_ + skinThickness_; |
486 |
} |
487 |
|
488 |
|
489 |
void InteractionManager::doPrePair(int *atid1, int *atid2, RealType *rij, RealType *rho_i_at_j, RealType *rho_j_at_i){ |
490 |
|
491 |
if (!initialized_) initialize(); |
492 |
|
493 |
DensityData ddat; |
494 |
|
495 |
ddat.atype1 = typeMap_[*atid1]; |
496 |
ddat.atype2 = typeMap_[*atid2]; |
497 |
ddat.rij = *rij; |
498 |
ddat.rho_i_at_j = *rho_i_at_j; |
499 |
ddat.rho_j_at_i = *rho_j_at_i; |
500 |
|
501 |
pair<AtomType*, AtomType*> key = make_pair(ddat.atype1, ddat.atype2); |
502 |
set<NonBondedInteraction*>::iterator it; |
503 |
|
504 |
for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){ |
505 |
if ((*it)->getFamily() == METALLIC_FAMILY) { |
506 |
dynamic_cast<MetallicInteraction*>(*it)->calcDensity(ddat); |
507 |
} |
508 |
} |
509 |
|
510 |
return; |
511 |
} |
512 |
|
513 |
void InteractionManager::doPreForce(int *atid, RealType *rho, RealType *frho, RealType *dfrhodrho){ |
514 |
|
515 |
if (!initialized_) initialize(); |
516 |
|
517 |
FunctionalData fdat; |
518 |
|
519 |
fdat.atype = typeMap_[*atid]; |
520 |
fdat.rho = *rho; |
521 |
fdat.frho = *frho; |
522 |
fdat.dfrhodrho = *dfrhodrho; |
523 |
|
524 |
pair<AtomType*, AtomType*> key = make_pair(fdat.atype, fdat.atype); |
525 |
set<NonBondedInteraction*>::iterator it; |
526 |
|
527 |
for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){ |
528 |
if ((*it)->getFamily() == METALLIC_FAMILY) { |
529 |
dynamic_cast<MetallicInteraction*>(*it)->calcFunctional(fdat); |
530 |
} |
531 |
} |
532 |
|
533 |
return; |
534 |
} |
535 |
|
536 |
void InteractionManager::doPair(int *atid1, int *atid2, RealType *d, RealType *r, RealType *r2, RealType *sw, int *topoDist, RealType *A1, RealType *A2, RealType *eFrame1, RealType *eFrame2, RealType *vpair, RealType *pot, RealType *f1, RealType *t1, RealType *t2, RealType *rho1, RealType *rho2, RealType *dfrho1, RealType *dfrho2, RealType *fshift1, RealType *fshift2){ |
537 |
|
538 |
if (!initialized_) initialize(); |
539 |
|
540 |
InteractionData idat; |
541 |
|
542 |
idat.atype1 = typeMap_[*atid1]; |
543 |
idat.atype2 = typeMap_[*atid2]; |
544 |
idat.d = Vector3d(d); |
545 |
idat.rij = *r; |
546 |
idat.r2 = *r2; |
547 |
idat.sw = *sw; |
548 |
idat.vdwMult = vdwScale_[*topoDist]; |
549 |
idat.electroMult = electrostaticScale_[*topoDist]; |
550 |
for (int i = 0; i < 4; i++) { |
551 |
idat.vpair[i] = vpair[i]; |
552 |
idat.pot[i] = pot[i]; |
553 |
} |
554 |
idat.f1 = Vector3d(f1[0], f1[1], f1[2]); |
555 |
idat.eFrame1 = Mat3x3d(eFrame1); |
556 |
idat.eFrame2 = Mat3x3d(eFrame2); |
557 |
idat.A1 = RotMat3x3d(A1); |
558 |
idat.A2 = RotMat3x3d(A2); |
559 |
idat.t1 = Vector3d(t1); |
560 |
idat.t2 = Vector3d(t2); |
561 |
idat.rho1 = *rho1; |
562 |
idat.rho2 = *rho2; |
563 |
idat.dfrho1 = *dfrho1; |
564 |
idat.dfrho2 = *dfrho2; |
565 |
idat.fshift1 = *fshift1; |
566 |
idat.fshift2 = *fshift2; |
567 |
|
568 |
pair<AtomType*, AtomType*> key = make_pair(idat.atype1, idat.atype2); |
569 |
set<NonBondedInteraction*>::iterator it; |
570 |
|
571 |
for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it) |
572 |
(*it)->calcForce(idat); |
573 |
|
574 |
for (int i = 0; i < 4; i++) { |
575 |
vpair[i] = idat.vpair[i]; |
576 |
pot[i] = idat.pot[i]; |
577 |
} |
578 |
|
579 |
for (int i = 0; i < 3; i++) { |
580 |
f1[i] = idat.f1[i]; |
581 |
t1[i] = idat.t1[i]; |
582 |
t2[i] = idat.t2[i]; |
583 |
} |
584 |
|
585 |
return; |
586 |
} |
587 |
|
588 |
void InteractionManager::doSkipCorrection(int *atid1, int *atid2, RealType *d, RealType *r, RealType *skippedCharge1, RealType *skippedCharge2, RealType *sw, RealType *electroMult, RealType *pot, RealType *vpair, RealType *f1, RealType *eFrame1, RealType *eFrame2, RealType *t1, RealType *t2){ |
589 |
|
590 |
if (!initialized_) initialize(); |
591 |
|
592 |
SkipCorrectionData skdat; |
593 |
|
594 |
skdat.atype1 = typeMap_[*atid1]; |
595 |
skdat.atype2 = typeMap_[*atid2]; |
596 |
skdat.d = Vector3d(d); |
597 |
skdat.rij = *r; |
598 |
skdat.skippedCharge1 = *skippedCharge1; |
599 |
skdat.skippedCharge2 = *skippedCharge2; |
600 |
skdat.sw = *sw; |
601 |
skdat.electroMult = *electroMult; |
602 |
for (int i = 0; i < 4; i++) { |
603 |
skdat.vpair[i] = vpair[i]; |
604 |
skdat.pot[i] = pot[i]; |
605 |
} |
606 |
skdat.f1 = Vector3d(f1); |
607 |
skdat.eFrame1 = Mat3x3d(eFrame1); |
608 |
skdat.eFrame2 = Mat3x3d(eFrame2); |
609 |
skdat.t1 = Vector3d(t1); |
610 |
skdat.t2 = Vector3d(t2); |
611 |
|
612 |
pair<AtomType*, AtomType*> key = make_pair(skdat.atype1, skdat.atype2); |
613 |
set<NonBondedInteraction*>::iterator it; |
614 |
|
615 |
for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){ |
616 |
if ((*it)->getFamily() == ELECTROSTATIC_FAMILY) { |
617 |
dynamic_cast<ElectrostaticInteraction*>(*it)->calcSkipCorrection(skdat); |
618 |
} |
619 |
} |
620 |
|
621 |
for (int i = 0; i < 4; i++) { |
622 |
vpair[i] = skdat.vpair[i]; |
623 |
pot[i] = skdat.pot[i]; |
624 |
} |
625 |
|
626 |
for (int i = 0; i < 3; i++) { |
627 |
f1[i] = skdat.f1[i]; |
628 |
t1[i] = skdat.t1[i]; |
629 |
t2[i] = skdat.t2[i]; |
630 |
} |
631 |
|
632 |
return; |
633 |
} |
634 |
|
635 |
void InteractionManager::doSelfCorrection(int *atid, RealType *eFrame, RealType *skippedCharge, RealType *pot, RealType *t){ |
636 |
|
637 |
if (!initialized_) initialize(); |
638 |
|
639 |
SelfCorrectionData scdat; |
640 |
|
641 |
scdat.atype = typeMap_[*atid]; |
642 |
scdat.eFrame = Mat3x3d(eFrame); |
643 |
scdat.skippedCharge = *skippedCharge; |
644 |
for (int i = 0; i < 4; i++){ |
645 |
scdat.pot[i] = pot[i]; |
646 |
} |
647 |
scdat.t = Vector3d(t); |
648 |
|
649 |
pair<AtomType*, AtomType*> key = make_pair(scdat.atype, scdat.atype); |
650 |
set<NonBondedInteraction*>::iterator it; |
651 |
|
652 |
for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it){ |
653 |
if ((*it)->getFamily() == ELECTROSTATIC_FAMILY) { |
654 |
dynamic_cast<ElectrostaticInteraction*>(*it)->calcSelfCorrection(scdat); |
655 |
} |
656 |
} |
657 |
|
658 |
t[0] = scdat.t.x(); |
659 |
t[1] = scdat.t.y(); |
660 |
t[2] = scdat.t.z(); |
661 |
|
662 |
return; |
663 |
} |
664 |
|
665 |
|
666 |
RealType InteractionManager::getSuggestedCutoffRadius(int *atid) { |
667 |
if (!initialized_) initialize(); |
668 |
|
669 |
AtomType* atype = typeMap_[*atid]; |
670 |
|
671 |
pair<AtomType*, AtomType*> key = make_pair(atype, atype); |
672 |
set<NonBondedInteraction*>::iterator it; |
673 |
RealType cutoff = 0.0; |
674 |
|
675 |
for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it) |
676 |
cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(atype, atype)); |
677 |
return cutoff; |
678 |
} |
679 |
|
680 |
RealType InteractionManager::getSuggestedCutoffRadius(AtomType* atype) { |
681 |
if (!initialized_) initialize(); |
682 |
|
683 |
pair<AtomType*, AtomType*> key = make_pair(atype, atype); |
684 |
set<NonBondedInteraction*>::iterator it; |
685 |
RealType cutoff = 0.0; |
686 |
|
687 |
for (it = interactions_[key].begin(); it != interactions_[key].end(); ++it) |
688 |
cutoff = max(cutoff, (*it)->getSuggestedCutoffRadius(atype, atype)); |
689 |
return cutoff; |
690 |
} |
691 |
|
692 |
void InteractionManager::getSwitch(RealType *r2, RealType *sw, RealType *dswdr, RealType *r, |
693 |
int *in_switching_region) { |
694 |
bool isr = switcher_->getSwitch( *r2 , *sw, *dswdr, *r); |
695 |
*in_switching_region = (int)isr; |
696 |
} |
697 |
|
698 |
} //end namespace OpenMD |
699 |
|
700 |
extern "C" { |
701 |
|
702 |
#define fortranDoPrePair FC_FUNC(do_prepair, DO_PREPAIR) |
703 |
#define fortranDoPreForce FC_FUNC(do_preforce, DO_PREFORCE) |
704 |
#define fortranDoPair FC_FUNC(do_pair, DO_PAIR) |
705 |
#define fortranDoSkipCorrection FC_FUNC(do_skip_correction, DO_SKIP_CORRECTION) |
706 |
#define fortranDoSelfCorrection FC_FUNC(do_self_correction, DO_SELF_CORRECTION) |
707 |
#define fortranGetCutoff FC_FUNC(get_cutoff, GET_CUTOFF) |
708 |
#define fortranGetSwitch FC_FUNC(get_switch, GET_SWITCH) |
709 |
|
710 |
void fortranDoPrePair(int *atid1, int *atid2, RealType *rij, |
711 |
RealType *rho_i_at_j, RealType *rho_j_at_i) { |
712 |
|
713 |
return OpenMD::InteractionManager::Instance()->doPrePair(atid1, atid2, rij, |
714 |
rho_i_at_j, |
715 |
rho_j_at_i); |
716 |
} |
717 |
void fortranDoPreForce(int *atid, RealType *rho, RealType *frho, |
718 |
RealType *dfrhodrho) { |
719 |
|
720 |
return OpenMD::InteractionManager::Instance()->doPreForce(atid, rho, frho, |
721 |
dfrhodrho); |
722 |
} |
723 |
|
724 |
void fortranDoPair(int *atid1, int *atid2, RealType *d, RealType *r, RealType *r2, |
725 |
RealType *sw, int *topoDist, RealType *A1, RealType *A2, |
726 |
RealType *eFrame1, RealType *eFrame2, |
727 |
RealType *vpair, RealType *pot, |
728 |
RealType *f1, RealType *t1, RealType *t2, |
729 |
RealType *rho1, RealType *rho2, RealType *dfrho1, RealType *dfrho2, |
730 |
RealType *fshift1, RealType *fshift2){ |
731 |
|
732 |
return OpenMD::InteractionManager::Instance()->doPair(atid1, atid2, d, r, r2, |
733 |
sw, topoDist, A1, A2, |
734 |
eFrame1, eFrame2, |
735 |
vpair, pot, |
736 |
f1, t1, t2, |
737 |
rho1, rho2, |
738 |
dfrho1, dfrho2, |
739 |
fshift1, fshift2); |
740 |
} |
741 |
|
742 |
void fortranDoSkipCorrection(int *atid1, int *atid2, RealType *d, |
743 |
RealType *r, RealType *skippedCharge1, |
744 |
RealType *skippedCharge2, RealType *sw, |
745 |
RealType *electroMult, RealType *pot, |
746 |
RealType *vpair, RealType *f1, |
747 |
RealType *eFrame1, RealType *eFrame2, |
748 |
RealType *t1, RealType *t2){ |
749 |
|
750 |
return OpenMD::InteractionManager::Instance()->doSkipCorrection(atid1, |
751 |
atid2, d, |
752 |
r, |
753 |
skippedCharge1, |
754 |
skippedCharge2, |
755 |
sw, electroMult, pot, |
756 |
vpair, f1, eFrame1, |
757 |
eFrame2, t1, t2); |
758 |
} |
759 |
|
760 |
void fortranDoSelfCorrection(int *atid, RealType *eFrame, RealType *skippedCharge, |
761 |
RealType *pot, RealType *t) { |
762 |
|
763 |
return OpenMD::InteractionManager::Instance()->doSelfCorrection(atid, |
764 |
eFrame, |
765 |
skippedCharge, |
766 |
pot, t); |
767 |
} |
768 |
RealType fortranGetCutoff(int *atid) { |
769 |
return OpenMD::InteractionManager::Instance()->getSuggestedCutoffRadius(atid); |
770 |
} |
771 |
|
772 |
void fortranGetSwitch(RealType *r2, RealType *sw, RealType *dswdr, RealType *r, |
773 |
int *in_switching_region) { |
774 |
return OpenMD::InteractionManager::Instance()->getSwitch(r2, sw, dswdr, r, |
775 |
in_switching_region); |
776 |
} |
777 |
} |