ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/nonbonded/Electrostatic.cpp
Revision: 1933
Committed: Fri Aug 23 15:59:23 2013 UTC (11 years, 8 months ago) by gezelter
File size: 54097 byte(s)
Log Message:
Fixed a bug (I think) that was caused by a cross product name collision.

File Contents

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

Properties

Name Value
svn:eol-style native