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 |
#ifdef IS_MPI |
44 |
#include <mpi.h> |
45 |
#endif |
46 |
|
47 |
#include <stdio.h> |
48 |
#include <string.h> |
49 |
|
50 |
#include <cmath> |
51 |
#include <numeric> |
52 |
#include "nonbonded/Electrostatic.hpp" |
53 |
#include "utils/simError.h" |
54 |
#include "types/NonBondedInteractionType.hpp" |
55 |
#include "types/FixedChargeAdapter.hpp" |
56 |
#include "types/FluctuatingChargeAdapter.hpp" |
57 |
#include "types/MultipoleAdapter.hpp" |
58 |
#include "io/Globals.hpp" |
59 |
#include "nonbonded/SlaterIntegrals.hpp" |
60 |
#include "utils/PhysicalConstants.hpp" |
61 |
#include "math/erfc.hpp" |
62 |
#include "math/SquareMatrix.hpp" |
63 |
#include "primitives/Molecule.hpp" |
64 |
#include "flucq/FluctuatingChargeForces.hpp" |
65 |
|
66 |
namespace OpenMD { |
67 |
|
68 |
Electrostatic::Electrostatic(): name_("Electrostatic"), initialized_(false), |
69 |
forceField_(NULL), info_(NULL), |
70 |
haveCutoffRadius_(false), |
71 |
haveDampingAlpha_(false), |
72 |
haveDielectric_(false), |
73 |
haveElectroSplines_(false) |
74 |
{ |
75 |
flucQ_ = new FluctuatingChargeForces(info_); |
76 |
} |
77 |
|
78 |
void Electrostatic::setForceField(ForceField *ff) { |
79 |
forceField_ = ff; |
80 |
flucQ_->setForceField(forceField_); |
81 |
} |
82 |
|
83 |
void Electrostatic::setSimulatedAtomTypes(set<AtomType*> &simtypes) { |
84 |
simTypes_ = simtypes; |
85 |
flucQ_->setSimulatedAtomTypes(simTypes_); |
86 |
} |
87 |
|
88 |
void Electrostatic::initialize() { |
89 |
|
90 |
Globals* simParams_ = info_->getSimParams(); |
91 |
|
92 |
summationMap_["HARD"] = esm_HARD; |
93 |
summationMap_["NONE"] = esm_HARD; |
94 |
summationMap_["SWITCHING_FUNCTION"] = esm_SWITCHING_FUNCTION; |
95 |
summationMap_["SHIFTED_POTENTIAL"] = esm_SHIFTED_POTENTIAL; |
96 |
summationMap_["SHIFTED_FORCE"] = esm_SHIFTED_FORCE; |
97 |
summationMap_["TAYLOR_SHIFTED"] = esm_TAYLOR_SHIFTED; |
98 |
summationMap_["REACTION_FIELD"] = esm_REACTION_FIELD; |
99 |
summationMap_["EWALD_FULL"] = esm_EWALD_FULL; |
100 |
summationMap_["EWALD_PME"] = esm_EWALD_PME; |
101 |
summationMap_["EWALD_SPME"] = esm_EWALD_SPME; |
102 |
screeningMap_["DAMPED"] = DAMPED; |
103 |
screeningMap_["UNDAMPED"] = UNDAMPED; |
104 |
|
105 |
// these prefactors convert the multipole interactions into kcal / mol |
106 |
// all were computed assuming distances are measured in angstroms |
107 |
// Charge-Charge, assuming charges are measured in electrons |
108 |
pre11_ = 332.0637778; |
109 |
// Charge-Dipole, assuming charges are measured in electrons, and |
110 |
// dipoles are measured in debyes |
111 |
pre12_ = 69.13373; |
112 |
// Dipole-Dipole, assuming dipoles are measured in Debye |
113 |
pre22_ = 14.39325; |
114 |
// Charge-Quadrupole, assuming charges are measured in electrons, and |
115 |
// quadrupoles are measured in 10^-26 esu cm^2 |
116 |
// This unit is also known affectionately as an esu centibarn. |
117 |
pre14_ = 69.13373; |
118 |
// Dipole-Quadrupole, assuming dipoles are measured in debyes and |
119 |
// quadrupoles in esu centibarns: |
120 |
pre24_ = 14.39325; |
121 |
// Quadrupole-Quadrupole, assuming esu centibarns: |
122 |
pre44_ = 14.39325; |
123 |
|
124 |
// conversions for the simulation box dipole moment |
125 |
chargeToC_ = 1.60217733e-19; |
126 |
angstromToM_ = 1.0e-10; |
127 |
debyeToCm_ = 3.33564095198e-30; |
128 |
|
129 |
// Default number of points for electrostatic splines |
130 |
np_ = 100; |
131 |
|
132 |
// variables to handle different summation methods for long-range |
133 |
// electrostatics: |
134 |
summationMethod_ = esm_HARD; |
135 |
screeningMethod_ = UNDAMPED; |
136 |
dielectric_ = 1.0; |
137 |
|
138 |
// check the summation method: |
139 |
if (simParams_->haveElectrostaticSummationMethod()) { |
140 |
string myMethod = simParams_->getElectrostaticSummationMethod(); |
141 |
toUpper(myMethod); |
142 |
map<string, ElectrostaticSummationMethod>::iterator i; |
143 |
i = summationMap_.find(myMethod); |
144 |
if ( i != summationMap_.end() ) { |
145 |
summationMethod_ = (*i).second; |
146 |
} else { |
147 |
// throw error |
148 |
sprintf( painCave.errMsg, |
149 |
"Electrostatic::initialize: Unknown electrostaticSummationMethod.\n" |
150 |
"\t(Input file specified %s .)\n" |
151 |
"\telectrostaticSummationMethod must be one of: \"hard\",\n" |
152 |
"\t\"shifted_potential\", \"shifted_force\",\n" |
153 |
"\t\"taylor_shifted\", or \"reaction_field\".\n", |
154 |
myMethod.c_str() ); |
155 |
painCave.isFatal = 1; |
156 |
simError(); |
157 |
} |
158 |
} else { |
159 |
// set ElectrostaticSummationMethod to the cutoffMethod: |
160 |
if (simParams_->haveCutoffMethod()){ |
161 |
string myMethod = simParams_->getCutoffMethod(); |
162 |
toUpper(myMethod); |
163 |
map<string, ElectrostaticSummationMethod>::iterator i; |
164 |
i = summationMap_.find(myMethod); |
165 |
if ( i != summationMap_.end() ) { |
166 |
summationMethod_ = (*i).second; |
167 |
} |
168 |
} |
169 |
} |
170 |
|
171 |
if (summationMethod_ == esm_REACTION_FIELD) { |
172 |
if (!simParams_->haveDielectric()) { |
173 |
// throw warning |
174 |
sprintf( painCave.errMsg, |
175 |
"SimInfo warning: dielectric was not specified in the input file\n\tfor the reaction field correction method.\n" |
176 |
"\tA default value of %f will be used for the dielectric.\n", dielectric_); |
177 |
painCave.isFatal = 0; |
178 |
painCave.severity = OPENMD_INFO; |
179 |
simError(); |
180 |
} else { |
181 |
dielectric_ = simParams_->getDielectric(); |
182 |
} |
183 |
haveDielectric_ = true; |
184 |
} |
185 |
|
186 |
if (simParams_->haveElectrostaticScreeningMethod()) { |
187 |
string myScreen = simParams_->getElectrostaticScreeningMethod(); |
188 |
toUpper(myScreen); |
189 |
map<string, ElectrostaticScreeningMethod>::iterator i; |
190 |
i = screeningMap_.find(myScreen); |
191 |
if ( i != screeningMap_.end()) { |
192 |
screeningMethod_ = (*i).second; |
193 |
} else { |
194 |
sprintf( painCave.errMsg, |
195 |
"SimInfo error: Unknown electrostaticScreeningMethod.\n" |
196 |
"\t(Input file specified %s .)\n" |
197 |
"\telectrostaticScreeningMethod must be one of: \"undamped\"\n" |
198 |
"or \"damped\".\n", myScreen.c_str() ); |
199 |
painCave.isFatal = 1; |
200 |
simError(); |
201 |
} |
202 |
} |
203 |
|
204 |
// check to make sure a cutoff value has been set: |
205 |
if (!haveCutoffRadius_) { |
206 |
sprintf( painCave.errMsg, "Electrostatic::initialize has no Default " |
207 |
"Cutoff value!\n"); |
208 |
painCave.severity = OPENMD_ERROR; |
209 |
painCave.isFatal = 1; |
210 |
simError(); |
211 |
} |
212 |
|
213 |
if (screeningMethod_ == DAMPED || summationMethod_ == esm_EWALD_FULL) { |
214 |
if (!simParams_->haveDampingAlpha()) { |
215 |
// first set a cutoff dependent alpha value |
216 |
// we assume alpha depends linearly with rcut from 0 to 20.5 ang |
217 |
dampingAlpha_ = 0.425 - cutoffRadius_* 0.02; |
218 |
if (dampingAlpha_ < 0.0) dampingAlpha_ = 0.0; |
219 |
// throw warning |
220 |
sprintf( painCave.errMsg, |
221 |
"Electrostatic::initialize: dampingAlpha was not specified in the\n" |
222 |
"\tinput file. A default value of %f (1/ang) will be used for the\n" |
223 |
"\tcutoff of %f (ang).\n", |
224 |
dampingAlpha_, cutoffRadius_); |
225 |
painCave.severity = OPENMD_INFO; |
226 |
painCave.isFatal = 0; |
227 |
simError(); |
228 |
} else { |
229 |
dampingAlpha_ = simParams_->getDampingAlpha(); |
230 |
} |
231 |
haveDampingAlpha_ = true; |
232 |
} |
233 |
|
234 |
|
235 |
Etypes.clear(); |
236 |
Etids.clear(); |
237 |
FQtypes.clear(); |
238 |
FQtids.clear(); |
239 |
ElectrostaticMap.clear(); |
240 |
Jij.clear(); |
241 |
nElectro_ = 0; |
242 |
nFlucq_ = 0; |
243 |
|
244 |
Etids.resize( forceField_->getNAtomType(), -1); |
245 |
FQtids.resize( forceField_->getNAtomType(), -1); |
246 |
|
247 |
set<AtomType*>::iterator at; |
248 |
for (at = simTypes_.begin(); at != simTypes_.end(); ++at) { |
249 |
if ((*at)->isElectrostatic()) nElectro_++; |
250 |
if ((*at)->isFluctuatingCharge()) nFlucq_++; |
251 |
} |
252 |
|
253 |
Jij.resize(nFlucq_); |
254 |
|
255 |
for (at = simTypes_.begin(); at != simTypes_.end(); ++at) { |
256 |
if ((*at)->isElectrostatic()) addType(*at); |
257 |
} |
258 |
|
259 |
if (summationMethod_ == esm_REACTION_FIELD) { |
260 |
preRF_ = (dielectric_ - 1.0) / |
261 |
((2.0 * dielectric_ + 1.0) * pow(cutoffRadius_,3) ); |
262 |
} |
263 |
|
264 |
RealType b0c, b1c, b2c, b3c, b4c, b5c; |
265 |
RealType db0c_1, db0c_2, db0c_3, db0c_4, db0c_5; |
266 |
RealType a2, expTerm, invArootPi; |
267 |
|
268 |
RealType r = cutoffRadius_; |
269 |
RealType r2 = r * r; |
270 |
RealType ric = 1.0 / r; |
271 |
RealType ric2 = ric * ric; |
272 |
|
273 |
if (screeningMethod_ == DAMPED) { |
274 |
a2 = dampingAlpha_ * dampingAlpha_; |
275 |
invArootPi = 1.0 / (dampingAlpha_ * sqrt(M_PI)); |
276 |
expTerm = exp(-a2 * r2); |
277 |
// values of Smith's B_l functions at the cutoff radius: |
278 |
b0c = erfc(dampingAlpha_ * r) / r; |
279 |
b1c = ( b0c + 2.0*a2 * expTerm * invArootPi) / r2; |
280 |
b2c = (3.0 * b1c + pow(2.0*a2, 2) * expTerm * invArootPi) / r2; |
281 |
b3c = (5.0 * b2c + pow(2.0*a2, 3) * expTerm * invArootPi) / r2; |
282 |
b4c = (7.0 * b3c + pow(2.0*a2, 4) * expTerm * invArootPi) / r2; |
283 |
b5c = (9.0 * b4c + pow(2.0*a2, 5) * expTerm * invArootPi) / r2; |
284 |
// Half the Smith self piece: |
285 |
selfMult1_ = - a2 * invArootPi; |
286 |
selfMult2_ = - 2.0 * a2 * a2 * invArootPi / 3.0; |
287 |
selfMult4_ = - 4.0 * a2 * a2 * a2 * invArootPi / 5.0; |
288 |
} else { |
289 |
a2 = 0.0; |
290 |
b0c = 1.0 / r; |
291 |
b1c = ( b0c) / r2; |
292 |
b2c = (3.0 * b1c) / r2; |
293 |
b3c = (5.0 * b2c) / r2; |
294 |
b4c = (7.0 * b3c) / r2; |
295 |
b5c = (9.0 * b4c) / r2; |
296 |
selfMult1_ = 0.0; |
297 |
selfMult2_ = 0.0; |
298 |
selfMult4_ = 0.0; |
299 |
} |
300 |
|
301 |
// higher derivatives of B_0 at the cutoff radius: |
302 |
db0c_1 = -r * b1c; |
303 |
db0c_2 = -b1c + r2 * b2c; |
304 |
db0c_3 = 3.0*r*b2c - r2*r*b3c; |
305 |
db0c_4 = 3.0*b2c - 6.0*r2*b3c + r2*r2*b4c; |
306 |
db0c_5 = -15.0*r*b3c + 10.0*r2*r*b4c - r2*r2*r*b5c; |
307 |
|
308 |
if (summationMethod_ != esm_EWALD_FULL) { |
309 |
selfMult1_ -= b0c; |
310 |
selfMult2_ += (db0c_2 + 2.0*db0c_1*ric) / 3.0; |
311 |
selfMult4_ -= (db0c_4 + 4.0*db0c_3*ric) / 15.0; |
312 |
} |
313 |
|
314 |
// working variables for the splines: |
315 |
RealType ri, ri2; |
316 |
RealType b0, b1, b2, b3, b4, b5; |
317 |
RealType db0_1, db0_2, db0_3, db0_4, db0_5; |
318 |
RealType f, fc, f0; |
319 |
RealType g, gc, g0, g1, g2, g3, g4; |
320 |
RealType h, hc, h1, h2, h3, h4; |
321 |
RealType s, sc, s2, s3, s4; |
322 |
RealType t, tc, t3, t4; |
323 |
RealType u, uc, u4; |
324 |
|
325 |
// working variables for Taylor expansion: |
326 |
RealType rmRc, rmRc2, rmRc3, rmRc4; |
327 |
|
328 |
// Approximate using splines using a maximum of 0.1 Angstroms |
329 |
// between points. |
330 |
int nptest = int((cutoffRadius_ + 2.0) / 0.1); |
331 |
np_ = (np_ > nptest) ? np_ : nptest; |
332 |
|
333 |
// Add a 2 angstrom safety window to deal with cutoffGroups that |
334 |
// have charged atoms longer than the cutoffRadius away from each |
335 |
// other. Splining is almost certainly the best choice here. |
336 |
// Direct calls to erfc would be preferrable if it is a very fast |
337 |
// implementation. |
338 |
|
339 |
RealType dx = (cutoffRadius_ + 2.0) / RealType(np_); |
340 |
|
341 |
// Storage vectors for the computed functions |
342 |
vector<RealType> rv; |
343 |
vector<RealType> v01v; |
344 |
vector<RealType> v11v; |
345 |
vector<RealType> v21v, v22v; |
346 |
vector<RealType> v31v, v32v; |
347 |
vector<RealType> v41v, v42v, v43v; |
348 |
|
349 |
for (int i = 1; i < np_ + 1; i++) { |
350 |
r = RealType(i) * dx; |
351 |
rv.push_back(r); |
352 |
|
353 |
ri = 1.0 / r; |
354 |
ri2 = ri * ri; |
355 |
|
356 |
r2 = r * r; |
357 |
expTerm = exp(-a2 * r2); |
358 |
|
359 |
// Taylor expansion factors (no need for factorials this way): |
360 |
rmRc = r - cutoffRadius_; |
361 |
rmRc2 = rmRc * rmRc / 2.0; |
362 |
rmRc3 = rmRc2 * rmRc / 3.0; |
363 |
rmRc4 = rmRc3 * rmRc / 4.0; |
364 |
|
365 |
// values of Smith's B_l functions at r: |
366 |
if (screeningMethod_ == DAMPED) { |
367 |
b0 = erfc(dampingAlpha_ * r) * ri; |
368 |
b1 = ( b0 + 2.0*a2 * expTerm * invArootPi) * ri2; |
369 |
b2 = (3.0 * b1 + pow(2.0*a2, 2) * expTerm * invArootPi) * ri2; |
370 |
b3 = (5.0 * b2 + pow(2.0*a2, 3) * expTerm * invArootPi) * ri2; |
371 |
b4 = (7.0 * b3 + pow(2.0*a2, 4) * expTerm * invArootPi) * ri2; |
372 |
b5 = (9.0 * b4 + pow(2.0*a2, 5) * expTerm * invArootPi) * ri2; |
373 |
} else { |
374 |
b0 = ri; |
375 |
b1 = ( b0) * ri2; |
376 |
b2 = (3.0 * b1) * ri2; |
377 |
b3 = (5.0 * b2) * ri2; |
378 |
b4 = (7.0 * b3) * ri2; |
379 |
b5 = (9.0 * b4) * ri2; |
380 |
} |
381 |
|
382 |
// higher derivatives of B_0 at r: |
383 |
db0_1 = -r * b1; |
384 |
db0_2 = -b1 + r2 * b2; |
385 |
db0_3 = 3.0*r*b2 - r2*r*b3; |
386 |
db0_4 = 3.0*b2 - 6.0*r2*b3 + r2*r2*b4; |
387 |
db0_5 = -15.0*r*b3 + 10.0*r2*r*b4 - r2*r2*r*b5; |
388 |
|
389 |
f = b0; |
390 |
fc = b0c; |
391 |
f0 = f - fc - rmRc*db0c_1; |
392 |
|
393 |
g = db0_1; |
394 |
gc = db0c_1; |
395 |
g0 = g - gc; |
396 |
g1 = g0 - rmRc *db0c_2; |
397 |
g2 = g1 - rmRc2*db0c_3; |
398 |
g3 = g2 - rmRc3*db0c_4; |
399 |
g4 = g3 - rmRc4*db0c_5; |
400 |
|
401 |
h = db0_2; |
402 |
hc = db0c_2; |
403 |
h1 = h - hc; |
404 |
h2 = h1 - rmRc *db0c_3; |
405 |
h3 = h2 - rmRc2*db0c_4; |
406 |
h4 = h3 - rmRc3*db0c_5; |
407 |
|
408 |
s = db0_3; |
409 |
sc = db0c_3; |
410 |
s2 = s - sc; |
411 |
s3 = s2 - rmRc *db0c_4; |
412 |
s4 = s3 - rmRc2*db0c_5; |
413 |
|
414 |
t = db0_4; |
415 |
tc = db0c_4; |
416 |
t3 = t - tc; |
417 |
t4 = t3 - rmRc *db0c_5; |
418 |
|
419 |
u = db0_5; |
420 |
uc = db0c_5; |
421 |
u4 = u - uc; |
422 |
|
423 |
// in what follows below, the various v functions are used for |
424 |
// potentials and torques, while the w functions show up in the |
425 |
// forces. |
426 |
|
427 |
switch (summationMethod_) { |
428 |
case esm_SHIFTED_FORCE: |
429 |
|
430 |
v01 = f - fc - rmRc*gc; |
431 |
v11 = g - gc - rmRc*hc; |
432 |
v21 = g*ri - gc*ric - rmRc*(hc - gc*ric)*ric; |
433 |
v22 = h - g*ri - (hc - gc*ric) - rmRc*(sc - (hc - gc*ric)*ric); |
434 |
v31 = (h-g*ri)*ri - (hc-gc*ric)*ric - rmRc*(sc-2.0*(hc-gc*ric)*ric)*ric; |
435 |
v32 = (s - 3.0*(h-g*ri)*ri) - (sc - 3.0*(hc-gc*ric)*ric) |
436 |
- rmRc*(tc - 3.0*(sc-2.0*(hc-gc*ric)*ric)*ric); |
437 |
v41 = (h - g*ri)*ri2 - (hc - gc*ric)*ric2 |
438 |
- rmRc*(sc - 3.0*(hc-gc*ric)*ric)*ric2; |
439 |
v42 = (s-3.0*(h-g*ri)*ri)*ri - (sc-3.0*(hc-gc*ric)*ric)*ric |
440 |
- rmRc*(tc - (4.0*sc - 9.0*(hc - gc*ric)*ric)*ric)*ric; |
441 |
|
442 |
v43 = (t - 3.0*(2.0*s - 5.0*(h - g*ri)*ri)*ri) |
443 |
- (tc - 3.0*(2.0*sc - 5.0*(hc - gc*ric)*ric)*ric) |
444 |
- rmRc*(uc-3.0*(2.0*tc - (7.0*sc - 15.0*(hc - gc*ric)*ric)*ric)*ric); |
445 |
|
446 |
dv01 = g - gc; |
447 |
dv11 = h - hc; |
448 |
dv21 = (h - g*ri)*ri - (hc - gc*ric)*ric; |
449 |
dv22 = (s - (h - g*ri)*ri) - (sc - (hc - gc*ric)*ric); |
450 |
dv31 = (s - 2.0*(h-g*ri)*ri)*ri - (sc - 2.0*(hc-gc*ric)*ric)*ric; |
451 |
dv32 = (t - 3.0*(s-2.0*(h-g*ri)*ri)*ri) |
452 |
- (tc - 3.0*(sc-2.0*(hc-gc*ric)*ric)*ric); |
453 |
dv41 = (s - 3.0*(h - g*ri)*ri)*ri2 - (sc - 3.0*(hc - gc*ric)*ric)*ric2; |
454 |
dv42 = (t - (4.0*s - 9.0*(h-g*ri)*ri)*ri)*ri |
455 |
- (tc - (4.0*sc - 9.0*(hc-gc*ric)*ric)*ric)*ric; |
456 |
dv43 = (u - 3.0*(2.0*t - (7.0*s - 15.0*(h - g*ri)*ri)*ri)*ri) |
457 |
- (uc - 3.0*(2.0*tc - (7.0*sc - 15.0*(hc - gc*ric)*ric)*ric)*ric); |
458 |
|
459 |
break; |
460 |
|
461 |
case esm_TAYLOR_SHIFTED: |
462 |
|
463 |
v01 = f0; |
464 |
v11 = g1; |
465 |
v21 = g2 * ri; |
466 |
v22 = h2 - v21; |
467 |
v31 = (h3 - g3 * ri) * ri; |
468 |
v32 = s3 - 3.0*v31; |
469 |
v41 = (h4 - g4 * ri) * ri2; |
470 |
v42 = s4 * ri - 3.0*v41; |
471 |
v43 = t4 - 6.0*v42 - 3.0*v41; |
472 |
|
473 |
dv01 = g0; |
474 |
dv11 = h1; |
475 |
dv21 = (h2 - g2*ri)*ri; |
476 |
dv22 = (s2 - (h2 - g2*ri)*ri); |
477 |
dv31 = (s3 - 2.0*(h3-g3*ri)*ri)*ri; |
478 |
dv32 = (t3 - 3.0*(s3-2.0*(h3-g3*ri)*ri)*ri); |
479 |
dv41 = (s4 - 3.0*(h4 - g4*ri)*ri)*ri2; |
480 |
dv42 = (t4 - (4.0*s4 - 9.0*(h4-g4*ri)*ri)*ri)*ri; |
481 |
dv43 = (u4 - 3.0*(2.0*t4 - (7.0*s4 - 15.0*(h4 - g4*ri)*ri)*ri)*ri); |
482 |
|
483 |
break; |
484 |
|
485 |
case esm_SHIFTED_POTENTIAL: |
486 |
|
487 |
v01 = f - fc; |
488 |
v11 = g - gc; |
489 |
v21 = g*ri - gc*ric; |
490 |
v22 = h - g*ri - (hc - gc*ric); |
491 |
v31 = (h-g*ri)*ri - (hc-gc*ric)*ric; |
492 |
v32 = (s - 3.0*(h-g*ri)*ri) - (sc - 3.0*(hc-gc*ric)*ric); |
493 |
v41 = (h - g*ri)*ri2 - (hc - gc*ric)*ric2; |
494 |
v42 = (s-3.0*(h-g*ri)*ri)*ri - (sc-3.0*(hc-gc*ric)*ric)*ric; |
495 |
v43 = (t - 3.0*(2.0*s - 5.0*(h - g*ri)*ri)*ri) |
496 |
- (tc - 3.0*(2.0*sc - 5.0*(hc - gc*ric)*ric)*ric); |
497 |
|
498 |
dv01 = g; |
499 |
dv11 = h; |
500 |
dv21 = (h - g*ri)*ri; |
501 |
dv22 = (s - (h - g*ri)*ri); |
502 |
dv31 = (s - 2.0*(h-g*ri)*ri)*ri; |
503 |
dv32 = (t - 3.0*(s-2.0*(h-g*ri)*ri)*ri); |
504 |
dv41 = (s - 3.0*(h - g*ri)*ri)*ri2; |
505 |
dv42 = (t - (4.0*s - 9.0*(h-g*ri)*ri)*ri)*ri; |
506 |
dv43 = (u - 3.0*(2.0*t - (7.0*s - 15.0*(h - g*ri)*ri)*ri)*ri); |
507 |
|
508 |
break; |
509 |
|
510 |
case esm_SWITCHING_FUNCTION: |
511 |
case esm_HARD: |
512 |
case esm_EWALD_FULL: |
513 |
|
514 |
v01 = f; |
515 |
v11 = g; |
516 |
v21 = g*ri; |
517 |
v22 = h - g*ri; |
518 |
v31 = (h-g*ri)*ri; |
519 |
v32 = (s - 3.0*(h-g*ri)*ri); |
520 |
v41 = (h - g*ri)*ri2; |
521 |
v42 = (s-3.0*(h-g*ri)*ri)*ri; |
522 |
v43 = (t - 3.0*(2.0*s - 5.0*(h - g*ri)*ri)*ri); |
523 |
|
524 |
dv01 = g; |
525 |
dv11 = h; |
526 |
dv21 = (h - g*ri)*ri; |
527 |
dv22 = (s - (h - g*ri)*ri); |
528 |
dv31 = (s - 2.0*(h-g*ri)*ri)*ri; |
529 |
dv32 = (t - 3.0*(s-2.0*(h-g*ri)*ri)*ri); |
530 |
dv41 = (s - 3.0*(h - g*ri)*ri)*ri2; |
531 |
dv42 = (t - (4.0*s - 9.0*(h-g*ri)*ri)*ri)*ri; |
532 |
dv43 = (u - 3.0*(2.0*t - (7.0*s - 15.0*(h - g*ri)*ri)*ri)*ri); |
533 |
|
534 |
break; |
535 |
|
536 |
case esm_REACTION_FIELD: |
537 |
|
538 |
// following DL_POLY's lead for shifting the image charge potential: |
539 |
f = b0 + preRF_ * r2; |
540 |
fc = b0c + preRF_ * cutoffRadius_ * cutoffRadius_; |
541 |
|
542 |
g = db0_1 + preRF_ * 2.0 * r; |
543 |
gc = db0c_1 + preRF_ * 2.0 * cutoffRadius_; |
544 |
|
545 |
h = db0_2 + preRF_ * 2.0; |
546 |
hc = db0c_2 + preRF_ * 2.0; |
547 |
|
548 |
v01 = f - fc; |
549 |
v11 = g - gc; |
550 |
v21 = g*ri - gc*ric; |
551 |
v22 = h - g*ri - (hc - gc*ric); |
552 |
v31 = (h-g*ri)*ri - (hc-gc*ric)*ric; |
553 |
v32 = (s - 3.0*(h-g*ri)*ri) - (sc - 3.0*(hc-gc*ric)*ric); |
554 |
v41 = (h - g*ri)*ri2 - (hc - gc*ric)*ric2; |
555 |
v42 = (s-3.0*(h-g*ri)*ri)*ri - (sc-3.0*(hc-gc*ric)*ric)*ric; |
556 |
v43 = (t - 3.0*(2.0*s - 5.0*(h - g*ri)*ri)*ri) |
557 |
- (tc - 3.0*(2.0*sc - 5.0*(hc - gc*ric)*ric)*ric); |
558 |
|
559 |
dv01 = g; |
560 |
dv11 = h; |
561 |
dv21 = (h - g*ri)*ri; |
562 |
dv22 = (s - (h - g*ri)*ri); |
563 |
dv31 = (s - 2.0*(h-g*ri)*ri)*ri; |
564 |
dv32 = (t - 3.0*(s-2.0*(h-g*ri)*ri)*ri); |
565 |
dv41 = (s - 3.0*(h - g*ri)*ri)*ri2; |
566 |
dv42 = (t - (4.0*s - 9.0*(h-g*ri)*ri)*ri)*ri; |
567 |
dv43 = (u - 3.0*(2.0*t - (7.0*s - 15.0*(h - g*ri)*ri)*ri)*ri); |
568 |
|
569 |
break; |
570 |
|
571 |
case esm_EWALD_PME: |
572 |
case esm_EWALD_SPME: |
573 |
default : |
574 |
map<string, ElectrostaticSummationMethod>::iterator i; |
575 |
std::string meth; |
576 |
for (i = summationMap_.begin(); i != summationMap_.end(); ++i) { |
577 |
if ((*i).second == summationMethod_) meth = (*i).first; |
578 |
} |
579 |
sprintf( painCave.errMsg, |
580 |
"Electrostatic::initialize: electrostaticSummationMethod %s \n" |
581 |
"\thas not been implemented yet. Please select one of:\n" |
582 |
"\t\"hard\", \"shifted_potential\", or \"shifted_force\"\n", |
583 |
meth.c_str() ); |
584 |
painCave.isFatal = 1; |
585 |
simError(); |
586 |
break; |
587 |
} |
588 |
|
589 |
// Add these computed values to the storage vectors for spline creation: |
590 |
v01v.push_back(v01); |
591 |
v11v.push_back(v11); |
592 |
v21v.push_back(v21); |
593 |
v22v.push_back(v22); |
594 |
v31v.push_back(v31); |
595 |
v32v.push_back(v32); |
596 |
v41v.push_back(v41); |
597 |
v42v.push_back(v42); |
598 |
v43v.push_back(v43); |
599 |
} |
600 |
|
601 |
// construct the spline structures and fill them with the values we've |
602 |
// computed: |
603 |
|
604 |
v01s = new CubicSpline(); |
605 |
v01s->addPoints(rv, v01v); |
606 |
v11s = new CubicSpline(); |
607 |
v11s->addPoints(rv, v11v); |
608 |
v21s = new CubicSpline(); |
609 |
v21s->addPoints(rv, v21v); |
610 |
v22s = new CubicSpline(); |
611 |
v22s->addPoints(rv, v22v); |
612 |
v31s = new CubicSpline(); |
613 |
v31s->addPoints(rv, v31v); |
614 |
v32s = new CubicSpline(); |
615 |
v32s->addPoints(rv, v32v); |
616 |
v41s = new CubicSpline(); |
617 |
v41s->addPoints(rv, v41v); |
618 |
v42s = new CubicSpline(); |
619 |
v42s->addPoints(rv, v42v); |
620 |
v43s = new CubicSpline(); |
621 |
v43s->addPoints(rv, v43v); |
622 |
|
623 |
haveElectroSplines_ = true; |
624 |
|
625 |
initialized_ = true; |
626 |
} |
627 |
|
628 |
void Electrostatic::addType(AtomType* atomType){ |
629 |
|
630 |
ElectrostaticAtomData electrostaticAtomData; |
631 |
electrostaticAtomData.is_Charge = false; |
632 |
electrostaticAtomData.is_Dipole = false; |
633 |
electrostaticAtomData.is_Quadrupole = false; |
634 |
electrostaticAtomData.is_Fluctuating = false; |
635 |
|
636 |
FixedChargeAdapter fca = FixedChargeAdapter(atomType); |
637 |
|
638 |
if (fca.isFixedCharge()) { |
639 |
electrostaticAtomData.is_Charge = true; |
640 |
electrostaticAtomData.fixedCharge = fca.getCharge(); |
641 |
} |
642 |
|
643 |
MultipoleAdapter ma = MultipoleAdapter(atomType); |
644 |
if (ma.isMultipole()) { |
645 |
if (ma.isDipole()) { |
646 |
electrostaticAtomData.is_Dipole = true; |
647 |
electrostaticAtomData.dipole = ma.getDipole(); |
648 |
} |
649 |
if (ma.isQuadrupole()) { |
650 |
electrostaticAtomData.is_Quadrupole = true; |
651 |
electrostaticAtomData.quadrupole = ma.getQuadrupole(); |
652 |
} |
653 |
} |
654 |
|
655 |
FluctuatingChargeAdapter fqa = FluctuatingChargeAdapter(atomType); |
656 |
|
657 |
if (fqa.isFluctuatingCharge()) { |
658 |
electrostaticAtomData.is_Fluctuating = true; |
659 |
electrostaticAtomData.electronegativity = fqa.getElectronegativity(); |
660 |
electrostaticAtomData.hardness = fqa.getHardness(); |
661 |
electrostaticAtomData.slaterN = fqa.getSlaterN(); |
662 |
electrostaticAtomData.slaterZeta = fqa.getSlaterZeta(); |
663 |
} |
664 |
|
665 |
int atid = atomType->getIdent(); |
666 |
int etid = Etypes.size(); |
667 |
int fqtid = FQtypes.size(); |
668 |
|
669 |
pair<set<int>::iterator,bool> ret; |
670 |
ret = Etypes.insert( atid ); |
671 |
if (ret.second == false) { |
672 |
sprintf( painCave.errMsg, |
673 |
"Electrostatic already had a previous entry with ident %d\n", |
674 |
atid); |
675 |
painCave.severity = OPENMD_INFO; |
676 |
painCave.isFatal = 0; |
677 |
simError(); |
678 |
} |
679 |
|
680 |
Etids[ atid ] = etid; |
681 |
ElectrostaticMap.push_back(electrostaticAtomData); |
682 |
|
683 |
if (electrostaticAtomData.is_Fluctuating) { |
684 |
ret = FQtypes.insert( atid ); |
685 |
if (ret.second == false) { |
686 |
sprintf( painCave.errMsg, |
687 |
"Electrostatic already had a previous fluctuating charge entry with ident %d\n", |
688 |
atid ); |
689 |
painCave.severity = OPENMD_INFO; |
690 |
painCave.isFatal = 0; |
691 |
simError(); |
692 |
} |
693 |
FQtids[atid] = fqtid; |
694 |
Jij[fqtid].resize(nFlucq_); |
695 |
|
696 |
// Now, iterate over all known fluctuating and add to the |
697 |
// coulomb integral map: |
698 |
|
699 |
std::set<int>::iterator it; |
700 |
for( it = FQtypes.begin(); it != FQtypes.end(); ++it) { |
701 |
int etid2 = Etids[ (*it) ]; |
702 |
int fqtid2 = FQtids[ (*it) ]; |
703 |
ElectrostaticAtomData eaData2 = ElectrostaticMap[ etid2 ]; |
704 |
RealType a = electrostaticAtomData.slaterZeta; |
705 |
RealType b = eaData2.slaterZeta; |
706 |
int m = electrostaticAtomData.slaterN; |
707 |
int n = eaData2.slaterN; |
708 |
|
709 |
// Create the spline of the coulombic integral for s-type |
710 |
// Slater orbitals. Add a 2 angstrom safety window to deal |
711 |
// with cutoffGroups that have charged atoms longer than the |
712 |
// cutoffRadius away from each other. |
713 |
|
714 |
RealType rval; |
715 |
RealType dr = (cutoffRadius_ + 2.0) / RealType(np_ - 1); |
716 |
vector<RealType> rvals; |
717 |
vector<RealType> Jvals; |
718 |
// don't start at i = 0, as rval = 0 is undefined for the |
719 |
// slater overlap integrals. |
720 |
for (int i = 1; i < np_+1; i++) { |
721 |
rval = RealType(i) * dr; |
722 |
rvals.push_back(rval); |
723 |
Jvals.push_back(sSTOCoulInt( a, b, m, n, rval * |
724 |
PhysicalConstants::angstromToBohr ) * |
725 |
PhysicalConstants::hartreeToKcal ); |
726 |
} |
727 |
|
728 |
CubicSpline* J = new CubicSpline(); |
729 |
J->addPoints(rvals, Jvals); |
730 |
Jij[fqtid][fqtid2] = J; |
731 |
Jij[fqtid2].resize( nFlucq_ ); |
732 |
Jij[fqtid2][fqtid] = J; |
733 |
} |
734 |
} |
735 |
return; |
736 |
} |
737 |
|
738 |
void Electrostatic::setCutoffRadius( RealType rCut ) { |
739 |
cutoffRadius_ = rCut; |
740 |
haveCutoffRadius_ = true; |
741 |
} |
742 |
|
743 |
void Electrostatic::setElectrostaticSummationMethod( ElectrostaticSummationMethod esm ) { |
744 |
summationMethod_ = esm; |
745 |
} |
746 |
void Electrostatic::setElectrostaticScreeningMethod( ElectrostaticScreeningMethod sm ) { |
747 |
screeningMethod_ = sm; |
748 |
} |
749 |
void Electrostatic::setDampingAlpha( RealType alpha ) { |
750 |
dampingAlpha_ = alpha; |
751 |
haveDampingAlpha_ = true; |
752 |
} |
753 |
void Electrostatic::setReactionFieldDielectric( RealType dielectric ){ |
754 |
dielectric_ = dielectric; |
755 |
haveDielectric_ = true; |
756 |
} |
757 |
|
758 |
void Electrostatic::calcForce(InteractionData &idat) { |
759 |
|
760 |
if (!initialized_) initialize(); |
761 |
|
762 |
data1 = ElectrostaticMap[Etids[idat.atid1]]; |
763 |
data2 = ElectrostaticMap[Etids[idat.atid2]]; |
764 |
|
765 |
U = 0.0; // Potential |
766 |
F.zero(); // Force |
767 |
Ta.zero(); // Torque on site a |
768 |
Tb.zero(); // Torque on site b |
769 |
Ea.zero(); // Electric field at site a |
770 |
Eb.zero(); // Electric field at site b |
771 |
Pa = 0.0; // Site potential at site a |
772 |
Pb = 0.0; // Site potential at site b |
773 |
dUdCa = 0.0; // fluctuating charge force at site a |
774 |
dUdCb = 0.0; // fluctuating charge force at site a |
775 |
|
776 |
// Indirect interactions mediated by the reaction field. |
777 |
indirect_Pot = 0.0; // Potential |
778 |
indirect_F.zero(); // Force |
779 |
indirect_Ta.zero(); // Torque on site a |
780 |
indirect_Tb.zero(); // Torque on site b |
781 |
|
782 |
// Excluded potential that is still computed for fluctuating charges |
783 |
excluded_Pot= 0.0; |
784 |
|
785 |
// some variables we'll need independent of electrostatic type: |
786 |
|
787 |
ri = 1.0 / *(idat.rij); |
788 |
rhat = *(idat.d) * ri; |
789 |
|
790 |
// logicals |
791 |
|
792 |
a_is_Charge = data1.is_Charge; |
793 |
a_is_Dipole = data1.is_Dipole; |
794 |
a_is_Quadrupole = data1.is_Quadrupole; |
795 |
a_is_Fluctuating = data1.is_Fluctuating; |
796 |
|
797 |
b_is_Charge = data2.is_Charge; |
798 |
b_is_Dipole = data2.is_Dipole; |
799 |
b_is_Quadrupole = data2.is_Quadrupole; |
800 |
b_is_Fluctuating = data2.is_Fluctuating; |
801 |
|
802 |
// Obtain all of the required radial function values from the |
803 |
// spline structures: |
804 |
|
805 |
// needed for fields (and forces): |
806 |
if (a_is_Charge || b_is_Charge) { |
807 |
v01s->getValueAndDerivativeAt( *(idat.rij), v01, dv01); |
808 |
} |
809 |
if (a_is_Dipole || b_is_Dipole) { |
810 |
v11s->getValueAndDerivativeAt( *(idat.rij), v11, dv11); |
811 |
v11or = ri * v11; |
812 |
} |
813 |
if (a_is_Quadrupole || b_is_Quadrupole || (a_is_Dipole && b_is_Dipole)) { |
814 |
v21s->getValueAndDerivativeAt( *(idat.rij), v21, dv21); |
815 |
v22s->getValueAndDerivativeAt( *(idat.rij), v22, dv22); |
816 |
v22or = ri * v22; |
817 |
} |
818 |
|
819 |
// needed for potentials (and forces and torques): |
820 |
if ((a_is_Dipole && b_is_Quadrupole) || |
821 |
(b_is_Dipole && a_is_Quadrupole)) { |
822 |
v31s->getValueAndDerivativeAt( *(idat.rij), v31, dv31); |
823 |
v32s->getValueAndDerivativeAt( *(idat.rij), v32, dv32); |
824 |
v31or = v31 * ri; |
825 |
v32or = v32 * ri; |
826 |
} |
827 |
if (a_is_Quadrupole && b_is_Quadrupole) { |
828 |
v41s->getValueAndDerivativeAt( *(idat.rij), v41, dv41); |
829 |
v42s->getValueAndDerivativeAt( *(idat.rij), v42, dv42); |
830 |
v43s->getValueAndDerivativeAt( *(idat.rij), v43, dv43); |
831 |
v42or = v42 * ri; |
832 |
v43or = v43 * ri; |
833 |
} |
834 |
|
835 |
// calculate the single-site contributions (fields, etc). |
836 |
|
837 |
if (a_is_Charge) { |
838 |
C_a = data1.fixedCharge; |
839 |
|
840 |
if (a_is_Fluctuating) { |
841 |
C_a += *(idat.flucQ1); |
842 |
} |
843 |
|
844 |
if (idat.excluded) { |
845 |
*(idat.skippedCharge2) += C_a; |
846 |
} else { |
847 |
// only do the field and site potentials if we're not excluded: |
848 |
Eb -= C_a * pre11_ * dv01 * rhat; |
849 |
Pb += C_a * pre11_ * v01; |
850 |
} |
851 |
} |
852 |
|
853 |
if (a_is_Dipole) { |
854 |
D_a = *(idat.dipole1); |
855 |
rdDa = dot(rhat, D_a); |
856 |
rxDa = cross(rhat, D_a); |
857 |
if (!idat.excluded) { |
858 |
Eb -= pre12_ * ((dv11-v11or) * rdDa * rhat + v11or * D_a); |
859 |
Pb += pre12_ * v11 * rdDa; |
860 |
} |
861 |
|
862 |
} |
863 |
|
864 |
if (a_is_Quadrupole) { |
865 |
Q_a = *(idat.quadrupole1); |
866 |
trQa = Q_a.trace(); |
867 |
Qar = Q_a * rhat; |
868 |
rQa = rhat * Q_a; |
869 |
rdQar = dot(rhat, Qar); |
870 |
rxQar = cross(rhat, Qar); |
871 |
if (!idat.excluded) { |
872 |
Eb -= pre14_ * (trQa * rhat * dv21 + 2.0 * Qar * v22or |
873 |
+ rdQar * rhat * (dv22 - 2.0*v22or)); |
874 |
Pb += pre14_ * (v21 * trQa + v22 * rdQar); |
875 |
} |
876 |
} |
877 |
|
878 |
if (b_is_Charge) { |
879 |
C_b = data2.fixedCharge; |
880 |
|
881 |
if (b_is_Fluctuating) |
882 |
C_b += *(idat.flucQ2); |
883 |
|
884 |
if (idat.excluded) { |
885 |
*(idat.skippedCharge1) += C_b; |
886 |
} else { |
887 |
// only do the field if we're not excluded: |
888 |
Ea += C_b * pre11_ * dv01 * rhat; |
889 |
Pa += C_b * pre11_ * v01; |
890 |
|
891 |
} |
892 |
} |
893 |
|
894 |
if (b_is_Dipole) { |
895 |
D_b = *(idat.dipole2); |
896 |
rdDb = dot(rhat, D_b); |
897 |
rxDb = cross(rhat, D_b); |
898 |
if (!idat.excluded) { |
899 |
Ea += pre12_ * ((dv11-v11or) * rdDb * rhat + v11or * D_b); |
900 |
Pa += pre12_ * v11 * rdDb; |
901 |
} |
902 |
} |
903 |
|
904 |
if (b_is_Quadrupole) { |
905 |
Q_b = *(idat.quadrupole2); |
906 |
trQb = Q_b.trace(); |
907 |
Qbr = Q_b * rhat; |
908 |
rQb = rhat * Q_b; |
909 |
rdQbr = dot(rhat, Qbr); |
910 |
rxQbr = cross(rhat, Qbr); |
911 |
if (!idat.excluded) { |
912 |
Ea += pre14_ * (trQb * rhat * dv21 + 2.0 * Qbr * v22or |
913 |
+ rdQbr * rhat * (dv22 - 2.0*v22or)); |
914 |
Pa += pre14_ * (v21 * trQb + v22 * rdQbr); |
915 |
} |
916 |
} |
917 |
|
918 |
|
919 |
if ((a_is_Fluctuating || b_is_Fluctuating) && idat.excluded) { |
920 |
J = Jij[FQtids[idat.atid1]][FQtids[idat.atid2]]; |
921 |
} |
922 |
|
923 |
if (a_is_Charge) { |
924 |
|
925 |
if (b_is_Charge) { |
926 |
pref = pre11_ * *(idat.electroMult); |
927 |
U += C_a * C_b * pref * v01; |
928 |
F += C_a * C_b * pref * dv01 * rhat; |
929 |
|
930 |
// If this is an excluded pair, there are still indirect |
931 |
// interactions via the reaction field we must worry about: |
932 |
|
933 |
if (summationMethod_ == esm_REACTION_FIELD && idat.excluded) { |
934 |
rfContrib = preRF_ * pref * C_a * C_b * *(idat.r2); |
935 |
indirect_Pot += rfContrib; |
936 |
indirect_F += rfContrib * 2.0 * ri * rhat; |
937 |
} |
938 |
|
939 |
// Fluctuating charge forces are handled via Coulomb integrals |
940 |
// for excluded pairs (i.e. those connected via bonds) and |
941 |
// with the standard charge-charge interaction otherwise. |
942 |
|
943 |
if (idat.excluded) { |
944 |
if (a_is_Fluctuating || b_is_Fluctuating) { |
945 |
coulInt = J->getValueAt( *(idat.rij) ); |
946 |
if (a_is_Fluctuating) dUdCa += C_b * coulInt; |
947 |
if (b_is_Fluctuating) dUdCb += C_a * coulInt; |
948 |
} |
949 |
} else { |
950 |
if (a_is_Fluctuating) dUdCa += C_b * pref * v01; |
951 |
if (b_is_Fluctuating) dUdCb += C_a * pref * v01; |
952 |
} |
953 |
} |
954 |
|
955 |
if (b_is_Dipole) { |
956 |
pref = pre12_ * *(idat.electroMult); |
957 |
U += C_a * pref * v11 * rdDb; |
958 |
F += C_a * pref * ((dv11 - v11or) * rdDb * rhat + v11or * D_b); |
959 |
Tb += C_a * pref * v11 * rxDb; |
960 |
|
961 |
if (a_is_Fluctuating) dUdCa += pref * v11 * rdDb; |
962 |
|
963 |
// Even if we excluded this pair from direct interactions, we |
964 |
// still have the reaction-field-mediated charge-dipole |
965 |
// interaction: |
966 |
|
967 |
if (summationMethod_ == esm_REACTION_FIELD && idat.excluded) { |
968 |
rfContrib = C_a * pref * preRF_ * 2.0 * *(idat.rij); |
969 |
indirect_Pot += rfContrib * rdDb; |
970 |
indirect_F += rfContrib * D_b / (*idat.rij); |
971 |
indirect_Tb += C_a * pref * preRF_ * rxDb; |
972 |
} |
973 |
} |
974 |
|
975 |
if (b_is_Quadrupole) { |
976 |
pref = pre14_ * *(idat.electroMult); |
977 |
U += C_a * pref * (v21 * trQb + v22 * rdQbr); |
978 |
F += C_a * pref * (trQb * dv21 * rhat + 2.0 * Qbr * v22or); |
979 |
F += C_a * pref * rdQbr * rhat * (dv22 - 2.0*v22or); |
980 |
Tb += C_a * pref * 2.0 * rxQbr * v22; |
981 |
|
982 |
if (a_is_Fluctuating) dUdCa += pref * (v21 * trQb + v22 * rdQbr); |
983 |
} |
984 |
} |
985 |
|
986 |
if (a_is_Dipole) { |
987 |
|
988 |
if (b_is_Charge) { |
989 |
pref = pre12_ * *(idat.electroMult); |
990 |
|
991 |
U -= C_b * pref * v11 * rdDa; |
992 |
F -= C_b * pref * ((dv11-v11or) * rdDa * rhat + v11or * D_a); |
993 |
Ta -= C_b * pref * v11 * rxDa; |
994 |
|
995 |
if (b_is_Fluctuating) dUdCb -= pref * v11 * rdDa; |
996 |
|
997 |
// Even if we excluded this pair from direct interactions, |
998 |
// we still have the reaction-field-mediated charge-dipole |
999 |
// interaction: |
1000 |
if (summationMethod_ == esm_REACTION_FIELD && idat.excluded) { |
1001 |
rfContrib = C_b * pref * preRF_ * 2.0 * *(idat.rij); |
1002 |
indirect_Pot -= rfContrib * rdDa; |
1003 |
indirect_F -= rfContrib * D_a / (*idat.rij); |
1004 |
indirect_Ta -= C_b * pref * preRF_ * rxDa; |
1005 |
} |
1006 |
} |
1007 |
|
1008 |
if (b_is_Dipole) { |
1009 |
pref = pre22_ * *(idat.electroMult); |
1010 |
DadDb = dot(D_a, D_b); |
1011 |
DaxDb = cross(D_a, D_b); |
1012 |
|
1013 |
U -= pref * (DadDb * v21 + rdDa * rdDb * v22); |
1014 |
F -= pref * (dv21 * DadDb * rhat + v22or * (rdDb * D_a + rdDa * D_b)); |
1015 |
F -= pref * (rdDa * rdDb) * (dv22 - 2.0*v22or) * rhat; |
1016 |
Ta += pref * ( v21 * DaxDb - v22 * rdDb * rxDa); |
1017 |
Tb += pref * (-v21 * DaxDb - v22 * rdDa * rxDb); |
1018 |
// Even if we excluded this pair from direct interactions, we |
1019 |
// still have the reaction-field-mediated dipole-dipole |
1020 |
// interaction: |
1021 |
if (summationMethod_ == esm_REACTION_FIELD && idat.excluded) { |
1022 |
rfContrib = -pref * preRF_ * 2.0; |
1023 |
indirect_Pot += rfContrib * DadDb; |
1024 |
indirect_Ta += rfContrib * DaxDb; |
1025 |
indirect_Tb -= rfContrib * DaxDb; |
1026 |
} |
1027 |
} |
1028 |
|
1029 |
if (b_is_Quadrupole) { |
1030 |
pref = pre24_ * *(idat.electroMult); |
1031 |
DadQb = D_a * Q_b; |
1032 |
DadQbr = dot(D_a, Qbr); |
1033 |
DaxQbr = cross(D_a, Qbr); |
1034 |
|
1035 |
U -= pref * ((trQb*rdDa + 2.0*DadQbr)*v31 + rdDa*rdQbr*v32); |
1036 |
F -= pref * (trQb*D_a + 2.0*DadQb) * v31or; |
1037 |
F -= pref * (trQb*rdDa + 2.0*DadQbr) * (dv31-v31or) * rhat; |
1038 |
F -= pref * (D_a*rdQbr + 2.0*rdDa*rQb) * v32or; |
1039 |
F -= pref * (rdDa * rdQbr * rhat * (dv32-3.0*v32or)); |
1040 |
Ta += pref * ((-trQb*rxDa + 2.0 * DaxQbr)*v31 - rxDa*rdQbr*v32); |
1041 |
Tb += pref * ((2.0*cross(DadQb, rhat) - 2.0*DaxQbr)*v31 |
1042 |
- 2.0*rdDa*rxQbr*v32); |
1043 |
} |
1044 |
} |
1045 |
|
1046 |
if (a_is_Quadrupole) { |
1047 |
if (b_is_Charge) { |
1048 |
pref = pre14_ * *(idat.electroMult); |
1049 |
U += C_b * pref * (v21 * trQa + v22 * rdQar); |
1050 |
F += C_b * pref * (trQa * rhat * dv21 + 2.0 * Qar * v22or); |
1051 |
F += C_b * pref * rdQar * rhat * (dv22 - 2.0*v22or); |
1052 |
Ta += C_b * pref * 2.0 * rxQar * v22; |
1053 |
|
1054 |
if (b_is_Fluctuating) dUdCb += pref * (v21 * trQa + v22 * rdQar); |
1055 |
} |
1056 |
if (b_is_Dipole) { |
1057 |
pref = pre24_ * *(idat.electroMult); |
1058 |
DbdQa = D_b * Q_a; |
1059 |
DbdQar = dot(D_b, Qar); |
1060 |
DbxQar = cross(D_b, Qar); |
1061 |
|
1062 |
U += pref * ((trQa*rdDb + 2.0*DbdQar)*v31 + rdDb*rdQar*v32); |
1063 |
F += pref * (trQa*D_b + 2.0*DbdQa) * v31or; |
1064 |
F += pref * (trQa*rdDb + 2.0*DbdQar) * (dv31-v31or) * rhat; |
1065 |
F += pref * (D_b*rdQar + 2.0*rdDb*rQa) * v32or; |
1066 |
F += pref * (rdDb * rdQar * rhat * (dv32-3.0*v32or)); |
1067 |
Ta += pref * ((-2.0*cross(DbdQa, rhat) + 2.0*DbxQar)*v31 |
1068 |
+ 2.0*rdDb*rxQar*v32); |
1069 |
Tb += pref * ((trQa*rxDb - 2.0 * DbxQar)*v31 + rxDb*rdQar*v32); |
1070 |
} |
1071 |
if (b_is_Quadrupole) { |
1072 |
pref = pre44_ * *(idat.electroMult); // yes |
1073 |
QaQb = Q_a * Q_b; |
1074 |
trQaQb = QaQb.trace(); |
1075 |
rQaQb = rhat * QaQb; |
1076 |
QaQbr = QaQb * rhat; |
1077 |
QaxQb = mCross(Q_a, Q_b); |
1078 |
rQaQbr = dot(rQa, Qbr); |
1079 |
rQaxQbr = cross(rQa, Qbr); |
1080 |
|
1081 |
U += pref * (trQa * trQb + 2.0 * trQaQb) * v41; |
1082 |
U += pref * (trQa * rdQbr + trQb * rdQar + 4.0 * rQaQbr) * v42; |
1083 |
U += pref * (rdQar * rdQbr) * v43; |
1084 |
|
1085 |
F += pref * rhat * (trQa * trQb + 2.0 * trQaQb)*dv41; |
1086 |
F += pref*rhat*(trQa*rdQbr + trQb*rdQar + 4.0*rQaQbr)*(dv42-2.0*v42or); |
1087 |
F += pref * rhat * (rdQar * rdQbr)*(dv43 - 4.0*v43or); |
1088 |
|
1089 |
F += pref * 2.0 * (trQb*rQa + trQa*rQb) * v42or; |
1090 |
F += pref * 4.0 * (rQaQb + QaQbr) * v42or; |
1091 |
F += pref * 2.0 * (rQa*rdQbr + rdQar*rQb) * v43or; |
1092 |
|
1093 |
Ta += pref * (- 4.0 * QaxQb * v41); |
1094 |
Ta += pref * (- 2.0 * trQb * cross(rQa, rhat) |
1095 |
+ 4.0 * cross(rhat, QaQbr) |
1096 |
- 4.0 * rQaxQbr) * v42; |
1097 |
Ta += pref * 2.0 * cross(rhat,Qar) * rdQbr * v43; |
1098 |
|
1099 |
|
1100 |
Tb += pref * (+ 4.0 * QaxQb * v41); |
1101 |
Tb += pref * (- 2.0 * trQa * cross(rQb, rhat) |
1102 |
- 4.0 * cross(rQaQb, rhat) |
1103 |
+ 4.0 * rQaxQbr) * v42; |
1104 |
// Possible replacement for line 2 above: |
1105 |
// + 4.0 * cross(rhat, QbQar) |
1106 |
|
1107 |
Tb += pref * 2.0 * cross(rhat,Qbr) * rdQar * v43; |
1108 |
} |
1109 |
} |
1110 |
|
1111 |
if (idat.doElectricField) { |
1112 |
*(idat.eField1) += Ea * *(idat.electroMult); |
1113 |
*(idat.eField2) += Eb * *(idat.electroMult); |
1114 |
} |
1115 |
|
1116 |
if (idat.doSitePotential) { |
1117 |
*(idat.sPot1) += Pa * *(idat.electroMult); |
1118 |
*(idat.sPot2) += Pb * *(idat.electroMult); |
1119 |
} |
1120 |
|
1121 |
if (a_is_Fluctuating) *(idat.dVdFQ1) += dUdCa * *(idat.sw); |
1122 |
if (b_is_Fluctuating) *(idat.dVdFQ2) += dUdCb * *(idat.sw); |
1123 |
|
1124 |
if (!idat.excluded) { |
1125 |
|
1126 |
*(idat.vpair) += U; |
1127 |
(*(idat.pot))[ELECTROSTATIC_FAMILY] += U * *(idat.sw); |
1128 |
*(idat.f1) += F * *(idat.sw); |
1129 |
|
1130 |
if (a_is_Dipole || a_is_Quadrupole) |
1131 |
*(idat.t1) += Ta * *(idat.sw); |
1132 |
|
1133 |
if (b_is_Dipole || b_is_Quadrupole) |
1134 |
*(idat.t2) += Tb * *(idat.sw); |
1135 |
|
1136 |
} else { |
1137 |
|
1138 |
// only accumulate the forces and torques resulting from the |
1139 |
// indirect reaction field terms. |
1140 |
|
1141 |
*(idat.vpair) += indirect_Pot; |
1142 |
(*(idat.excludedPot))[ELECTROSTATIC_FAMILY] += excluded_Pot; |
1143 |
(*(idat.pot))[ELECTROSTATIC_FAMILY] += *(idat.sw) * indirect_Pot; |
1144 |
*(idat.f1) += *(idat.sw) * indirect_F; |
1145 |
|
1146 |
if (a_is_Dipole || a_is_Quadrupole) |
1147 |
*(idat.t1) += *(idat.sw) * indirect_Ta; |
1148 |
|
1149 |
if (b_is_Dipole || b_is_Quadrupole) |
1150 |
*(idat.t2) += *(idat.sw) * indirect_Tb; |
1151 |
} |
1152 |
return; |
1153 |
} |
1154 |
|
1155 |
void Electrostatic::calcSelfCorrection(SelfData &sdat) { |
1156 |
|
1157 |
if (!initialized_) initialize(); |
1158 |
|
1159 |
ElectrostaticAtomData data = ElectrostaticMap[Etids[sdat.atid]]; |
1160 |
|
1161 |
// logicals |
1162 |
bool i_is_Charge = data.is_Charge; |
1163 |
bool i_is_Dipole = data.is_Dipole; |
1164 |
bool i_is_Quadrupole = data.is_Quadrupole; |
1165 |
bool i_is_Fluctuating = data.is_Fluctuating; |
1166 |
RealType C_a = data.fixedCharge; |
1167 |
RealType self(0.0), preVal, DdD, trQ, trQQ; |
1168 |
|
1169 |
if (i_is_Dipole) { |
1170 |
DdD = data.dipole.lengthSquare(); |
1171 |
} |
1172 |
|
1173 |
if (i_is_Fluctuating) { |
1174 |
C_a += *(sdat.flucQ); |
1175 |
|
1176 |
flucQ_->getSelfInteraction(sdat.atid, *(sdat.flucQ), |
1177 |
(*(sdat.excludedPot))[ELECTROSTATIC_FAMILY], |
1178 |
*(sdat.flucQfrc) ); |
1179 |
|
1180 |
} |
1181 |
|
1182 |
switch (summationMethod_) { |
1183 |
case esm_REACTION_FIELD: |
1184 |
|
1185 |
if (i_is_Charge) { |
1186 |
// Self potential [see Wang and Hermans, "Reaction Field |
1187 |
// Molecular Dynamics Simulation with Friedman’s Image Charge |
1188 |
// Method," J. Phys. Chem. 99, 12001-12007 (1995).] |
1189 |
preVal = pre11_ * preRF_ * C_a * C_a; |
1190 |
(*(sdat.pot))[ELECTROSTATIC_FAMILY] -= 0.5 * preVal / cutoffRadius_; |
1191 |
} |
1192 |
|
1193 |
if (i_is_Dipole) { |
1194 |
(*(sdat.pot))[ELECTROSTATIC_FAMILY] -= pre22_ * preRF_ * DdD; |
1195 |
} |
1196 |
|
1197 |
break; |
1198 |
|
1199 |
case esm_SHIFTED_FORCE: |
1200 |
case esm_SHIFTED_POTENTIAL: |
1201 |
case esm_TAYLOR_SHIFTED: |
1202 |
case esm_EWALD_FULL: |
1203 |
if (i_is_Charge) |
1204 |
self += selfMult1_ * pre11_ * C_a * (C_a + *(sdat.skippedCharge)); |
1205 |
if (i_is_Dipole) |
1206 |
self += selfMult2_ * pre22_ * DdD; |
1207 |
if (i_is_Quadrupole) { |
1208 |
trQ = data.quadrupole.trace(); |
1209 |
trQQ = (data.quadrupole * data.quadrupole).trace(); |
1210 |
self += selfMult4_ * pre44_ * (2.0*trQQ + trQ*trQ); |
1211 |
if (i_is_Charge) |
1212 |
self -= selfMult2_ * pre14_ * 2.0 * C_a * trQ; |
1213 |
} |
1214 |
(*(sdat.pot))[ELECTROSTATIC_FAMILY] += self; |
1215 |
break; |
1216 |
default: |
1217 |
break; |
1218 |
} |
1219 |
} |
1220 |
|
1221 |
RealType Electrostatic::getSuggestedCutoffRadius(pair<AtomType*, AtomType*> atypes) { |
1222 |
// This seems to work moderately well as a default. There's no |
1223 |
// inherent scale for 1/r interactions that we can standardize. |
1224 |
// 12 angstroms seems to be a reasonably good guess for most |
1225 |
// cases. |
1226 |
return 12.0; |
1227 |
} |
1228 |
|
1229 |
|
1230 |
void Electrostatic::ReciprocalSpaceSum(RealType& pot) { |
1231 |
|
1232 |
RealType kPot = 0.0; |
1233 |
RealType kVir = 0.0; |
1234 |
|
1235 |
const RealType mPoleConverter = 0.20819434; // converts from the |
1236 |
// internal units of |
1237 |
// Debye (for dipoles) |
1238 |
// or Debye-angstroms |
1239 |
// (for quadrupoles) to |
1240 |
// electron angstroms or |
1241 |
// electron-angstroms^2 |
1242 |
|
1243 |
const RealType eConverter = 332.0637778; // convert the |
1244 |
// Charge-Charge |
1245 |
// electrostatic |
1246 |
// interactions into kcal / |
1247 |
// mol assuming distances |
1248 |
// are measured in |
1249 |
// angstroms. |
1250 |
|
1251 |
Mat3x3d hmat = info_->getSnapshotManager()->getCurrentSnapshot()->getHmat(); |
1252 |
Vector3d box = hmat.diagonals(); |
1253 |
RealType boxMax = box.max(); |
1254 |
|
1255 |
//int kMax = int(2.0 * M_PI / (pow(dampingAlpha_,2)*cutoffRadius_ * boxMax) ); |
1256 |
int kMax = 7; |
1257 |
int kSqMax = kMax*kMax + 2; |
1258 |
|
1259 |
int kLimit = kMax+1; |
1260 |
int kLim2 = 2*kMax+1; |
1261 |
int kSqLim = kSqMax; |
1262 |
|
1263 |
vector<RealType> AK(kSqLim+1, 0.0); |
1264 |
RealType xcl = 2.0 * M_PI / box.x(); |
1265 |
RealType ycl = 2.0 * M_PI / box.y(); |
1266 |
RealType zcl = 2.0 * M_PI / box.z(); |
1267 |
RealType rcl = 2.0 * M_PI / boxMax; |
1268 |
RealType rvol = 2.0 * M_PI /(box.x() * box.y() * box.z()); |
1269 |
|
1270 |
if(dampingAlpha_ < 1.0e-12) return; |
1271 |
|
1272 |
RealType ralph = -0.25/pow(dampingAlpha_,2); |
1273 |
|
1274 |
// Calculate and store exponential factors |
1275 |
|
1276 |
vector<vector<RealType> > elc; |
1277 |
vector<vector<RealType> > emc; |
1278 |
vector<vector<RealType> > enc; |
1279 |
vector<vector<RealType> > els; |
1280 |
vector<vector<RealType> > ems; |
1281 |
vector<vector<RealType> > ens; |
1282 |
|
1283 |
int nMax = info_->getNAtoms(); |
1284 |
|
1285 |
elc.resize(kLimit+1); |
1286 |
emc.resize(kLimit+1); |
1287 |
enc.resize(kLimit+1); |
1288 |
els.resize(kLimit+1); |
1289 |
ems.resize(kLimit+1); |
1290 |
ens.resize(kLimit+1); |
1291 |
|
1292 |
for (int j = 0; j < kLimit+1; j++) { |
1293 |
elc[j].resize(nMax); |
1294 |
emc[j].resize(nMax); |
1295 |
enc[j].resize(nMax); |
1296 |
els[j].resize(nMax); |
1297 |
ems[j].resize(nMax); |
1298 |
ens[j].resize(nMax); |
1299 |
} |
1300 |
|
1301 |
Vector3d t( 2.0 * M_PI ); |
1302 |
t.Vdiv(t, box); |
1303 |
|
1304 |
SimInfo::MoleculeIterator mi; |
1305 |
Molecule::AtomIterator ai; |
1306 |
int i; |
1307 |
Vector3d r; |
1308 |
Vector3d tt; |
1309 |
|
1310 |
for (Molecule* mol = info_->beginMolecule(mi); mol != NULL; |
1311 |
mol = info_->nextMolecule(mi)) { |
1312 |
for(Atom* atom = mol->beginAtom(ai); atom != NULL; |
1313 |
atom = mol->nextAtom(ai)) { |
1314 |
|
1315 |
i = atom->getLocalIndex(); |
1316 |
r = atom->getPos(); |
1317 |
info_->getSnapshotManager()->getCurrentSnapshot()->wrapVector(r); |
1318 |
|
1319 |
tt.Vmul(t, r); |
1320 |
|
1321 |
elc[1][i] = 1.0; |
1322 |
emc[1][i] = 1.0; |
1323 |
enc[1][i] = 1.0; |
1324 |
els[1][i] = 0.0; |
1325 |
ems[1][i] = 0.0; |
1326 |
ens[1][i] = 0.0; |
1327 |
|
1328 |
elc[2][i] = cos(tt.x()); |
1329 |
emc[2][i] = cos(tt.y()); |
1330 |
enc[2][i] = cos(tt.z()); |
1331 |
els[2][i] = sin(tt.x()); |
1332 |
ems[2][i] = sin(tt.y()); |
1333 |
ens[2][i] = sin(tt.z()); |
1334 |
|
1335 |
for(int l = 3; l <= kLimit; l++) { |
1336 |
elc[l][i]=elc[l-1][i]*elc[2][i]-els[l-1][i]*els[2][i]; |
1337 |
emc[l][i]=emc[l-1][i]*emc[2][i]-ems[l-1][i]*ems[2][i]; |
1338 |
enc[l][i]=enc[l-1][i]*enc[2][i]-ens[l-1][i]*ens[2][i]; |
1339 |
els[l][i]=els[l-1][i]*elc[2][i]+elc[l-1][i]*els[2][i]; |
1340 |
ems[l][i]=ems[l-1][i]*emc[2][i]+emc[l-1][i]*ems[2][i]; |
1341 |
ens[l][i]=ens[l-1][i]*enc[2][i]+enc[l-1][i]*ens[2][i]; |
1342 |
} |
1343 |
} |
1344 |
} |
1345 |
|
1346 |
// Calculate and store AK coefficients: |
1347 |
|
1348 |
RealType eksq = 1.0; |
1349 |
RealType expf = 0.0; |
1350 |
if (ralph < 0.0) expf = exp(ralph*rcl*rcl); |
1351 |
for (i = 1; i <= kSqLim; i++) { |
1352 |
RealType rksq = float(i)*rcl*rcl; |
1353 |
eksq = expf*eksq; |
1354 |
AK[i] = eConverter * eksq/rksq; |
1355 |
} |
1356 |
|
1357 |
/* |
1358 |
* Loop over all k vectors k = 2 pi (ll/Lx, mm/Ly, nn/Lz) |
1359 |
* the values of ll, mm and nn are selected so that the symmetry of |
1360 |
* reciprocal lattice is taken into account i.e. the following |
1361 |
* rules apply. |
1362 |
* |
1363 |
* ll ranges over the values 0 to kMax only. |
1364 |
* |
1365 |
* mm ranges over 0 to kMax when ll=0 and over |
1366 |
* -kMax to kMax otherwise. |
1367 |
* nn ranges over 1 to kMax when ll=mm=0 and over |
1368 |
* -kMax to kMax otherwise. |
1369 |
* |
1370 |
* Hence the result of the summation must be doubled at the end. |
1371 |
*/ |
1372 |
|
1373 |
std::vector<RealType> clm(nMax, 0.0); |
1374 |
std::vector<RealType> slm(nMax, 0.0); |
1375 |
std::vector<RealType> ckr(nMax, 0.0); |
1376 |
std::vector<RealType> skr(nMax, 0.0); |
1377 |
std::vector<RealType> ckc(nMax, 0.0); |
1378 |
std::vector<RealType> cks(nMax, 0.0); |
1379 |
std::vector<RealType> dkc(nMax, 0.0); |
1380 |
std::vector<RealType> dks(nMax, 0.0); |
1381 |
std::vector<RealType> qkc(nMax, 0.0); |
1382 |
std::vector<RealType> qks(nMax, 0.0); |
1383 |
std::vector<Vector3d> dxk(nMax, V3Zero); |
1384 |
std::vector<Vector3d> qxk(nMax, V3Zero); |
1385 |
RealType rl, rm, rn; |
1386 |
Vector3d kVec; |
1387 |
Vector3d Qk; |
1388 |
Mat3x3d k2; |
1389 |
RealType ckcs, ckss, dkcs, dkss, qkcs, qkss; |
1390 |
int atid; |
1391 |
ElectrostaticAtomData data; |
1392 |
RealType C, dk, qk; |
1393 |
Vector3d D; |
1394 |
Mat3x3d Q; |
1395 |
|
1396 |
int mMin = kLimit; |
1397 |
int nMin = kLimit + 1; |
1398 |
for (int l = 1; l <= kLimit; l++) { |
1399 |
int ll = l - 1; |
1400 |
rl = xcl * float(ll); |
1401 |
for (int mmm = mMin; mmm <= kLim2; mmm++) { |
1402 |
int mm = mmm - kLimit; |
1403 |
int m = abs(mm) + 1; |
1404 |
rm = ycl * float(mm); |
1405 |
// Set temporary products of exponential terms |
1406 |
for (Molecule* mol = info_->beginMolecule(mi); mol != NULL; |
1407 |
mol = info_->nextMolecule(mi)) { |
1408 |
for(Atom* atom = mol->beginAtom(ai); atom != NULL; |
1409 |
atom = mol->nextAtom(ai)) { |
1410 |
|
1411 |
i = atom->getLocalIndex(); |
1412 |
if(mm < 0) { |
1413 |
clm[i]=elc[l][i]*emc[m][i]+els[l][i]*ems[m][i]; |
1414 |
slm[i]=els[l][i]*emc[m][i]-ems[m][i]*elc[l][i]; |
1415 |
} else { |
1416 |
clm[i]=elc[l][i]*emc[m][i]-els[l][i]*ems[m][i]; |
1417 |
slm[i]=els[l][i]*emc[m][i]+ems[m][i]*elc[l][i]; |
1418 |
} |
1419 |
} |
1420 |
} |
1421 |
for (int nnn = nMin; nnn <= kLim2; nnn++) { |
1422 |
int nn = nnn - kLimit; |
1423 |
int n = abs(nn) + 1; |
1424 |
rn = zcl * float(nn); |
1425 |
// Test on magnitude of k vector: |
1426 |
int kk=ll*ll + mm*mm + nn*nn; |
1427 |
if(kk <= kSqLim) { |
1428 |
kVec = Vector3d(rl, rm, rn); |
1429 |
k2 = outProduct(kVec, kVec); |
1430 |
// Calculate exp(ikr) terms |
1431 |
for (Molecule* mol = info_->beginMolecule(mi); mol != NULL; |
1432 |
mol = info_->nextMolecule(mi)) { |
1433 |
for(Atom* atom = mol->beginAtom(ai); atom != NULL; |
1434 |
atom = mol->nextAtom(ai)) { |
1435 |
i = atom->getLocalIndex(); |
1436 |
|
1437 |
if (nn < 0) { |
1438 |
ckr[i]=clm[i]*enc[n][i]+slm[i]*ens[n][i]; |
1439 |
skr[i]=slm[i]*enc[n][i]-clm[i]*ens[n][i]; |
1440 |
|
1441 |
} else { |
1442 |
ckr[i]=clm[i]*enc[n][i]-slm[i]*ens[n][i]; |
1443 |
skr[i]=slm[i]*enc[n][i]+clm[i]*ens[n][i]; |
1444 |
} |
1445 |
} |
1446 |
} |
1447 |
|
1448 |
// Calculate scalar and vector products for each site: |
1449 |
|
1450 |
for (Molecule* mol = info_->beginMolecule(mi); mol != NULL; |
1451 |
mol = info_->nextMolecule(mi)) { |
1452 |
for(Atom* atom = mol->beginAtom(ai); atom != NULL; |
1453 |
atom = mol->nextAtom(ai)) { |
1454 |
i = atom->getLocalIndex(); |
1455 |
int atid = atom->getAtomType()->getIdent(); |
1456 |
data = ElectrostaticMap[Etids[atid]]; |
1457 |
|
1458 |
if (data.is_Charge) { |
1459 |
C = data.fixedCharge; |
1460 |
if (atom->isFluctuatingCharge()) C += atom->getFlucQPos(); |
1461 |
ckc[i] = C * ckr[i]; |
1462 |
cks[i] = C * skr[i]; |
1463 |
} |
1464 |
|
1465 |
if (data.is_Dipole) { |
1466 |
D = atom->getDipole() * mPoleConverter; |
1467 |
dk = dot(D, kVec); |
1468 |
dxk[i] = cross(D, kVec); |
1469 |
dkc[i] = dk * ckr[i]; |
1470 |
dks[i] = dk * skr[i]; |
1471 |
} |
1472 |
if (data.is_Quadrupole) { |
1473 |
Q = atom->getQuadrupole() * mPoleConverter; |
1474 |
Qk = Q * kVec; |
1475 |
qk = dot(kVec, Qk); |
1476 |
qxk[i] = -cross(kVec, Qk); |
1477 |
qkc[i] = qk * ckr[i]; |
1478 |
qks[i] = qk * skr[i]; |
1479 |
} |
1480 |
} |
1481 |
} |
1482 |
|
1483 |
// calculate vector sums |
1484 |
|
1485 |
ckcs = std::accumulate(ckc.begin(),ckc.end(),0.0); |
1486 |
ckss = std::accumulate(cks.begin(),cks.end(),0.0); |
1487 |
dkcs = std::accumulate(dkc.begin(),dkc.end(),0.0); |
1488 |
dkss = std::accumulate(dks.begin(),dks.end(),0.0); |
1489 |
qkcs = std::accumulate(qkc.begin(),qkc.end(),0.0); |
1490 |
qkss = std::accumulate(qks.begin(),qks.end(),0.0); |
1491 |
|
1492 |
#ifdef IS_MPI |
1493 |
MPI_Allreduce(MPI_IN_PLACE, &ckcs, 1, MPI_REALTYPE, |
1494 |
MPI_SUM, MPI_COMM_WORLD); |
1495 |
MPI_Allreduce(MPI_IN_PLACE, &ckss, 1, MPI_REALTYPE, |
1496 |
MPI_SUM, MPI_COMM_WORLD); |
1497 |
MPI_Allreduce(MPI_IN_PLACE, &dkcs, 1, MPI_REALTYPE, |
1498 |
MPI_SUM, MPI_COMM_WORLD); |
1499 |
MPI_Allreduce(MPI_IN_PLACE, &dkss, 1, MPI_REALTYPE, |
1500 |
MPI_SUM, MPI_COMM_WORLD); |
1501 |
MPI_Allreduce(MPI_IN_PLACE, &qkcs, 1, MPI_REALTYPE, |
1502 |
MPI_SUM, MPI_COMM_WORLD); |
1503 |
MPI_Allreduce(MPI_IN_PLACE, &qkss, 1, MPI_REALTYPE, |
1504 |
MPI_SUM, MPI_COMM_WORLD); |
1505 |
#endif |
1506 |
|
1507 |
// Accumulate potential energy and virial contribution: |
1508 |
|
1509 |
kPot += 2.0 * rvol * AK[kk]*((ckss+dkcs-qkss)*(ckss+dkcs-qkss) |
1510 |
+ (ckcs-dkss-qkcs)*(ckcs-dkss-qkcs)); |
1511 |
|
1512 |
kVir += 2.0 * rvol * AK[kk]*(ckcs*ckcs+ckss*ckss |
1513 |
+4.0*(ckss*dkcs-ckcs*dkss) |
1514 |
+3.0*(dkcs*dkcs+dkss*dkss) |
1515 |
-6.0*(ckss*qkss+ckcs*qkcs) |
1516 |
+8.0*(dkss*qkcs-dkcs*qkss) |
1517 |
+5.0*(qkss*qkss+qkcs*qkcs)); |
1518 |
|
1519 |
// Calculate force and torque for each site: |
1520 |
|
1521 |
for (Molecule* mol = info_->beginMolecule(mi); mol != NULL; |
1522 |
mol = info_->nextMolecule(mi)) { |
1523 |
for(Atom* atom = mol->beginAtom(ai); atom != NULL; |
1524 |
atom = mol->nextAtom(ai)) { |
1525 |
|
1526 |
i = atom->getLocalIndex(); |
1527 |
atid = atom->getAtomType()->getIdent(); |
1528 |
data = ElectrostaticMap[Etids[atid]]; |
1529 |
|
1530 |
RealType qfrc = AK[kk]*((cks[i]+dkc[i]-qks[i])*(ckcs-dkss-qkcs) |
1531 |
- (ckc[i]-dks[i]-qkc[i])*(ckss+dkcs-qkss)); |
1532 |
RealType qtrq1 = AK[kk]*(skr[i]*(ckcs-dkss-qkcs) |
1533 |
-ckr[i]*(ckss+dkcs-qkss)); |
1534 |
RealType qtrq2 = 2.0*AK[kk]*(ckr[i]*(ckcs-dkss-qkcs) |
1535 |
+skr[i]*(ckss+dkcs-qkss)); |
1536 |
|
1537 |
atom->addFrc( 4.0 * rvol * qfrc * kVec ); |
1538 |
|
1539 |
if (atom->isFluctuatingCharge()) { |
1540 |
atom->addFlucQFrc( - 2.0 * rvol * qtrq2 ); |
1541 |
} |
1542 |
|
1543 |
if (data.is_Dipole) { |
1544 |
atom->addTrq( 4.0 * rvol * qtrq1 * dxk[i] ); |
1545 |
} |
1546 |
if (data.is_Quadrupole) { |
1547 |
atom->addTrq( 4.0 * rvol * qtrq2 * qxk[i] ); |
1548 |
} |
1549 |
} |
1550 |
} |
1551 |
} |
1552 |
} |
1553 |
nMin = 1; |
1554 |
} |
1555 |
mMin = 1; |
1556 |
} |
1557 |
pot += kPot; |
1558 |
} |
1559 |
} |